Stratimikos Version of the Day
Loading...
Searching...
No Matches
Thyra_BelosLinearOpWithSolveFactory_def.hpp
1/*
2// @HEADER
3// ***********************************************************************
4//
5// Stratimikos: Thyra-based strategies for linear solvers
6// Copyright (2006) Sandia Corporation
7//
8// Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
9// license for use of this work by or on behalf of the U.S. Government.
10//
11// Redistribution and use in source and binary forms, with or without
12// modification, are permitted provided that the following conditions are
13// met:
14//
15// 1. Redistributions of source code must retain the above copyright
16// notice, this list of conditions and the following disclaimer.
17//
18// 2. Redistributions in binary form must reproduce the above copyright
19// notice, this list of conditions and the following disclaimer in the
20// documentation and/or other materials provided with the distribution.
21//
22// 3. Neither the name of the Corporation nor the names of the
23// contributors may be used to endorse or promote products derived from
24// this software without specific prior written permission.
25//
26// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37//
38// Questions? Contact Roscoe A. Bartlett (rabartl@sandia.gov)
39//
40// ***********************************************************************
41// @HEADER
42*/
43
44
45#ifndef THYRA_BELOS_LINEAR_OP_WITH_SOLVE_FACTORY_HPP
46#define THYRA_BELOS_LINEAR_OP_WITH_SOLVE_FACTORY_HPP
47
48
49#include "Thyra_BelosLinearOpWithSolveFactory_decl.hpp"
50#include "Thyra_BelosLinearOpWithSolve.hpp"
51#include "Thyra_ScaledAdjointLinearOpBase.hpp"
52
58#include "BelosGCRODRSolMgr.hpp"
59#include "BelosRCGSolMgr.hpp"
60#include "BelosMinresSolMgr.hpp"
61#include "BelosTFQMRSolMgr.hpp"
64#include "BelosThyraAdapter.hpp"
65
66#include "Thyra_BelosTpetrasSolverAdapter.hpp"
67
68#include "Teuchos_VerboseObjectParameterListHelpers.hpp"
69#include "Teuchos_StandardParameterEntryValidators.hpp"
70#include "Teuchos_ParameterList.hpp"
71#include "Teuchos_dyn_cast.hpp"
72#include "Teuchos_ValidatorXMLConverterDB.hpp"
73#include "Teuchos_StandardValidatorXMLConverters.hpp"
74
75
76namespace Thyra {
77
78
79// Parameter names for Parameter List
80
81template<class Scalar>
82const std::string BelosLinearOpWithSolveFactory<Scalar>::SolverType_name = "Solver Type";
83template<class Scalar>
84const std::string BelosLinearOpWithSolveFactory<Scalar>::SolverType_default = "Pseudo Block GMRES";
85template<class Scalar>
86const std::string BelosLinearOpWithSolveFactory<Scalar>::SolverTypes_name = "Solver Types";
87template<class Scalar>
88const std::string BelosLinearOpWithSolveFactory<Scalar>::BlockGMRES_name = "Block GMRES";
89template<class Scalar>
90const std::string BelosLinearOpWithSolveFactory<Scalar>::PseudoBlockGMRES_name = "Pseudo Block GMRES";
91template<class Scalar>
92const std::string BelosLinearOpWithSolveFactory<Scalar>::BlockCG_name = "Block CG";
93template<class Scalar>
94const std::string BelosLinearOpWithSolveFactory<Scalar>::PseudoBlockCG_name = "Pseudo Block CG";
95template<class Scalar>
96const std::string BelosLinearOpWithSolveFactory<Scalar>::PseudoBlockStochasticCG_name = "Pseudo Block Stochastic CG";
97template<class Scalar>
99template<class Scalar>
101template<class Scalar>
102const std::string BelosLinearOpWithSolveFactory<Scalar>::MINRES_name = "MINRES";
103template<class Scalar>
104const std::string BelosLinearOpWithSolveFactory<Scalar>::TFQMR_name = "TFQMR";
105template<class Scalar>
106const std::string BelosLinearOpWithSolveFactory<Scalar>::BiCGStab_name = "BiCGStab";
107template<class Scalar>
108const std::string BelosLinearOpWithSolveFactory<Scalar>::FixedPoint_name = "Fixed Point";
109template<class Scalar>
110const std::string BelosLinearOpWithSolveFactory<Scalar>::TpetraGmres_name = "TPETRA GMRES";
111template<class Scalar>
112const std::string BelosLinearOpWithSolveFactory<Scalar>::TpetraGmresPipeline_name = "TPETRA GMRES PIPELINE";
113template<class Scalar>
114const std::string BelosLinearOpWithSolveFactory<Scalar>::TpetraGmresSingleReduce_name = "TPETRA GMRES SINGLE REDUCE";
115template<class Scalar>
116const std::string BelosLinearOpWithSolveFactory<Scalar>::TpetraGmresSstep_name = "TPETRA GMRES S-STEP";
117
118template<class Scalar>
119const std::string BelosLinearOpWithSolveFactory<Scalar>::ConvergenceTestFrequency_name = "Convergence Test Frequency";
120
121namespace {
122const std::string LeftPreconditionerIfUnspecified_name = "Left Preconditioner If Unspecified";
123}
124
125// Constructors/initializers/accessors
126
127
128template<class Scalar>
130 :solverType_(SOLVER_TYPE_PSEUDO_BLOCK_GMRES),
131 convergenceTestFrequency_(1)
132{
133 updateThisValidParamList();
134}
135
136
137template<class Scalar>
139 const RCP<PreconditionerFactoryBase<Scalar> > &precFactory
140 )
141 :solverType_(SOLVER_TYPE_PSEUDO_BLOCK_GMRES)
142{
143 this->setPreconditionerFactory(precFactory, "");
144}
145
146
147// Overridden from LinearOpWithSolveFactoryBase
148
149
150template<class Scalar>
155
156
157template<class Scalar>
159 const RCP<PreconditionerFactoryBase<Scalar> > &precFactory,
160 const std::string &precFactoryName
161 )
162{
163 TEUCHOS_TEST_FOR_EXCEPT(!precFactory.get());
164 RCP<const Teuchos::ParameterList>
165 precFactoryValidPL = precFactory->getValidParameters();
166 const std::string _precFactoryName =
167 ( precFactoryName != ""
168 ? precFactoryName
169 : ( precFactoryValidPL.get() ? precFactoryValidPL->name() : "GENERIC PRECONDITIONER FACTORY" )
170 );
171 precFactory_ = precFactory;
172 precFactoryName_ = _precFactoryName;
173 updateThisValidParamList();
174}
175
176
177template<class Scalar>
178RCP<PreconditionerFactoryBase<Scalar> >
183
184
185template<class Scalar>
187 RCP<PreconditionerFactoryBase<Scalar> > *precFactory,
188 std::string *precFactoryName
189 )
190{
191 if(precFactory) *precFactory = precFactory_;
192 if(precFactoryName) *precFactoryName = precFactoryName_;
193 precFactory_ = Teuchos::null;
194 precFactoryName_ = "";
195 updateThisValidParamList();
196}
197
198
199template<class Scalar>
201 const LinearOpSourceBase<Scalar> &fwdOpSrc
202 ) const
203{
204 if(precFactory_.get())
205 return precFactory_->isCompatible(fwdOpSrc);
206 return true; // Without a preconditioner, we are compatible with all linear operators!
207}
208
209
210template<class Scalar>
211RCP<LinearOpWithSolveBase<Scalar> >
216
217
218template<class Scalar>
220 const RCP<const LinearOpSourceBase<Scalar> > &fwdOpSrc,
221 LinearOpWithSolveBase<Scalar> *Op,
222 const ESupportSolveUse supportSolveUse
223 ) const
224{
225 using Teuchos::null;
226 initializeOpImpl(fwdOpSrc,null,null,false,Op,supportSolveUse);
227}
228
229
230template<class Scalar>
232 const RCP<const LinearOpSourceBase<Scalar> > &fwdOpSrc,
233 LinearOpWithSolveBase<Scalar> *Op
234 ) const
235{
236 using Teuchos::null;
237 initializeOpImpl(fwdOpSrc,null,null,true,Op,SUPPORT_SOLVE_UNSPECIFIED);
238}
239
240
241template<class Scalar>
243 const EPreconditionerInputType precOpType
244 ) const
245{
246 if(precFactory_.get())
247 return true;
248 return (precOpType==PRECONDITIONER_INPUT_TYPE_AS_OPERATOR);
249}
250
251
252template<class Scalar>
254 const RCP<const LinearOpSourceBase<Scalar> > &fwdOpSrc,
255 const RCP<const PreconditionerBase<Scalar> > &prec,
256 LinearOpWithSolveBase<Scalar> *Op,
257 const ESupportSolveUse supportSolveUse
258 ) const
259{
260 using Teuchos::null;
261 initializeOpImpl(fwdOpSrc,null,prec,false,Op,supportSolveUse);
262}
263
264
265template<class Scalar>
267 const RCP<const LinearOpSourceBase<Scalar> > &fwdOpSrc,
268 const RCP<const LinearOpSourceBase<Scalar> > &approxFwdOpSrc,
269 LinearOpWithSolveBase<Scalar> *Op,
270 const ESupportSolveUse supportSolveUse
271 ) const
272{
273 using Teuchos::null;
274 initializeOpImpl(fwdOpSrc,approxFwdOpSrc,null,false,Op,supportSolveUse);
275}
276
277
278template<class Scalar>
280 LinearOpWithSolveBase<Scalar> *Op,
281 RCP<const LinearOpSourceBase<Scalar> > *fwdOpSrc,
282 RCP<const PreconditionerBase<Scalar> > *prec,
283 RCP<const LinearOpSourceBase<Scalar> > *approxFwdOpSrc,
284 ESupportSolveUse *supportSolveUse
285 ) const
286{
287#ifdef TEUCHOS_DEBUG
288 TEUCHOS_TEST_FOR_EXCEPT(Op==NULL);
289#endif
291 &belosOp = Teuchos::dyn_cast<BelosLinearOpWithSolve<Scalar> >(*Op);
292 RCP<const LinearOpSourceBase<Scalar> >
293 _fwdOpSrc = belosOp.extract_fwdOpSrc();
294 RCP<const PreconditionerBase<Scalar> >
295 _prec = ( belosOp.isExternalPrec() ? belosOp.extract_prec() : Teuchos::null );
296 // Note: above we only extract the preconditioner if it was passed in
297 // externally. Otherwise, we need to hold on to it so that we can reuse it
298 // in the next initialization.
299 RCP<const LinearOpSourceBase<Scalar> >
300 _approxFwdOpSrc = belosOp.extract_approxFwdOpSrc();
301 ESupportSolveUse
302 _supportSolveUse = belosOp.supportSolveUse();
303 if(fwdOpSrc) *fwdOpSrc = _fwdOpSrc;
304 if(prec) *prec = _prec;
305 if(approxFwdOpSrc) *approxFwdOpSrc = _approxFwdOpSrc;
306 if(supportSolveUse) *supportSolveUse = _supportSolveUse;
307}
308
309
310// Overridden from ParameterListAcceptor
311
312
313template<class Scalar>
315 RCP<Teuchos::ParameterList> const& paramList
316 )
317{
318 TEUCHOS_TEST_FOR_EXCEPT(paramList.get()==NULL);
319 paramList->validateParametersAndSetDefaults(*this->getValidParameters(), 1);
320 paramList_ = paramList;
321 solverType_ =
322 Teuchos::getIntegralValue<EBelosSolverType>(*paramList_, SolverType_name);
323 convergenceTestFrequency_ =
324 Teuchos::getParameter<int>(*paramList_, ConvergenceTestFrequency_name);
325 Teuchos::readVerboseObjectSublist(&*paramList_,this);
326}
327
328
329template<class Scalar>
330RCP<Teuchos::ParameterList>
335
336
337template<class Scalar>
338RCP<Teuchos::ParameterList>
340{
341 RCP<Teuchos::ParameterList> _paramList = paramList_;
342 paramList_ = Teuchos::null;
343 return _paramList;
344}
345
346
347template<class Scalar>
348RCP<const Teuchos::ParameterList>
350{
351 return paramList_;
352}
353
354
355template<class Scalar>
356RCP<const Teuchos::ParameterList>
358{
359 return thisValidParamList_;
360}
361
362
363// Public functions overridden from Teuchos::Describable
364
365
366template<class Scalar>
368{
369 std::ostringstream oss;
370 oss << Teuchos::Describable::description();
371 oss << "{";
372 if (nonnull(precFactory_)) {
373 oss << precFactory_->description();
374 }
375 oss << "}";
376 return oss.str();
377}
378
379
380// private
381
382
383template<class Scalar>
384RCP<const Teuchos::ParameterList>
385BelosLinearOpWithSolveFactory<Scalar>::generateAndGetValidParameters()
386{
387 using Teuchos::as;
388 using Teuchos::tuple;
389 using Teuchos::setStringToIntegralParameter;
390Teuchos::ValidatorXMLConverterDB::addConverter(
391 Teuchos::DummyObjectGetter<
392 Teuchos::StringToIntegralParameterEntryValidator<EBelosSolverType>
393 >::getDummyObject(),
394 Teuchos::DummyObjectGetter<Teuchos::StringToIntegralValidatorXMLConverter<
395 EBelosSolverType> >::getDummyObject());
396
397 typedef MultiVectorBase<Scalar> MV_t;
398 typedef LinearOpBase<Scalar> LO_t;
399 static RCP<Teuchos::ParameterList> validParamList;
400 if(validParamList.get()==NULL) {
401 validParamList = Teuchos::rcp(new Teuchos::ParameterList("BelosLinearOpWithSolveFactory"));
402 setStringToIntegralParameter<EBelosSolverType>(
403 SolverType_name, SolverType_default,
404 "Type of linear solver algorithm to use.",
405 tuple<std::string>(
406 "Block GMRES",
407 "Pseudo Block GMRES",
408 "Block CG",
409 "Pseudo Block CG",
410 "Pseudo Block Stochastic CG",
411 "GCRODR",
412 "RCG",
413 "MINRES",
414 "TFQMR",
415 "BiCGStab",
416 "Fixed Point",
417 "TPETRA GMRES",
418 "TPETRA GMRES PIPELINE",
419 "TPETRA GMRES SINGLE REDUCE",
420 "TPETRA GMRES S-STEP"
421 ),
422 tuple<std::string>(
423 "Block GMRES solver for nonsymmetric linear systems. It can also solve "
424 "single right-hand side systems, and can also perform Flexible GMRES "
425 "(where the preconditioner may change at every iteration, for example "
426 "for inner-outer iterations), by setting options in the \"Block GMRES\" "
427 "sublist.",
428
429 "GMRES solver for nonsymmetric linear systems, that performs single "
430 "right-hand side solves on multiple right-hand sides at once. It "
431 "exploits operator multivector multiplication in order to amortize "
432 "global communication costs. Individual linear systems are deflated "
433 "out as they are solved.",
434
435 "Block CG solver for symmetric (Hermitian in complex arithmetic) "
436 "positive definite linear systems. It can also solve single "
437 "right-hand-side systems.",
438
439 "CG solver that performs single right-hand side CG on multiple right-hand "
440 "sides at once. It exploits operator multivector multiplication in order "
441 "to amortize global communication costs. Individual linear systems are "
442 "deflated out as they are solved.",
443
444 "stochastic CG solver that performs single right-hand side CG on multiple right-hand "
445 "sides at once. It exploits operator multivector multiplication in order "
446 "to amortize global communication costs. Individual linear systems are "
447 "deflated out as they are solved. [EXPERIMENTAL]",
448
449 "Variant of GMRES that performs subspace recycling to accelerate "
450 "convergence for sequences of solves with related linear systems. "
451 "Individual linear systems are deflated out as they are solved. "
452 "The current implementation only supports real-valued Scalar types.",
453
454 "CG solver for symmetric (Hermitian in complex arithmetic) positive "
455 "definite linear systems, that performs subspace recycling to "
456 "accelerate convergence for sequences of related linear systems.",
457
458 "MINRES solver for symmetric indefinite linear systems. It performs "
459 "single-right-hand-side solves on multiple right-hand sides sequentially.",
460
461 "TFQMR (Transpose-Free QMR) solver for nonsymmetric linear systems.",
462
463 "BiCGStab solver for nonsymmetric linear systems.",
464
465 "Fixed point iteration",
466
467 "Native Tpetra implementation of GMRES",
468
469 "Native Tpetra implementation of pipeline GMRES",
470
471 "Native Tpetra implementation of single-reduce GMRES",
472
473 "Native Tpetra implementation of s-step GMRES"
474 ),
475 tuple<EBelosSolverType>(
476 SOLVER_TYPE_BLOCK_GMRES,
477 SOLVER_TYPE_PSEUDO_BLOCK_GMRES,
478 SOLVER_TYPE_BLOCK_CG,
479 SOLVER_TYPE_PSEUDO_BLOCK_CG,
480 SOLVER_TYPE_PSEUDO_BLOCK_STOCHASTIC_CG,
481 SOLVER_TYPE_GCRODR,
482 SOLVER_TYPE_RCG,
483 SOLVER_TYPE_MINRES,
484 SOLVER_TYPE_TFQMR,
485 SOLVER_TYPE_BICGSTAB,
486 SOLVER_TYPE_FIXEDPOINT,
487 SOLVER_TYPE_TPETRA_GMRES,
488 SOLVER_TYPE_TPETRA_GMRES_PIPELINE,
489 SOLVER_TYPE_TPETRA_GMRES_SINGLE_REDUCE,
490 SOLVER_TYPE_TPETRA_GMRES_SSTEP
491 ),
492 &*validParamList
493 );
494 validParamList->set(ConvergenceTestFrequency_name, as<int>(1),
495 "Number of linear solver iterations to skip between applying"
496 " user-defined convergence test.");
497 validParamList->set(
498 LeftPreconditionerIfUnspecified_name, false,
499 "If the preconditioner does not specify if it is left or right, and this\n"
500 "option is set to true, put the preconditioner on the left side.\n"
501 "Historically, preconditioning is on the right. Some solvers may not\n"
502 "support left preconditioning.");
503 Teuchos::ParameterList
504 &solverTypesSL = validParamList->sublist(SolverTypes_name);
505
506 const bool lapackSupportsScalar = Belos::Details::LapackSupportsScalar<Scalar>::value;
507 const bool scalarIsReal = !Teuchos::ScalarTraits<Scalar>::isComplex;
508
509 {
511 solverTypesSL.sublist(BlockGMRES_name).setParameters(
512 *mgr.getValidParameters()
513 );
514 }
515 {
517 solverTypesSL.sublist(PseudoBlockGMRES_name).setParameters(
518 *mgr.getValidParameters()
519 );
520 }
521 if (lapackSupportsScalar) {
522 Belos::BlockCGSolMgr<Scalar,MV_t,LO_t> mgr;
523 solverTypesSL.sublist(BlockCG_name).setParameters(
524 *mgr.getValidParameters()
525 );
526 }
527 if (lapackSupportsScalar) {
528 Belos::PseudoBlockCGSolMgr<Scalar,MV_t,LO_t> mgr;
529 solverTypesSL.sublist(PseudoBlockCG_name).setParameters(
530 *mgr.getValidParameters()
531 );
532 }
533 {
534 Belos::PseudoBlockStochasticCGSolMgr<Scalar,MV_t,LO_t> mgr;
535 solverTypesSL.sublist(PseudoBlockStochasticCG_name).setParameters(
536 *mgr.getValidParameters()
537 );
538 }
539 if (lapackSupportsScalar) {
540 Belos::GCRODRSolMgr<Scalar,MV_t,LO_t> mgr;
541 solverTypesSL.sublist(GCRODR_name).setParameters(
542 *mgr.getValidParameters()
543 );
544 }
545 if (lapackSupportsScalar && scalarIsReal) {
546 Belos::RCGSolMgr<Scalar,MV_t,LO_t> mgr;
547 solverTypesSL.sublist(RCG_name).setParameters(
548 *mgr.getValidParameters()
549 );
550 }
551 {
552 Belos::MinresSolMgr<Scalar,MV_t,LO_t> mgr;
553 solverTypesSL.sublist(MINRES_name).setParameters(
554 *mgr.getValidParameters()
555 );
556 }
557 {
558 Belos::TFQMRSolMgr<Scalar,MV_t,LO_t> mgr;
559 solverTypesSL.sublist(TFQMR_name).setParameters(
560 *mgr.getValidParameters()
561 );
562 }
563 {
564 Belos::BiCGStabSolMgr<Scalar,MV_t,LO_t> mgr;
565 solverTypesSL.sublist(BiCGStab_name).setParameters(
566 *mgr.getValidParameters()
567 );
568 }
569 {
570 Belos::FixedPointSolMgr<Scalar,MV_t,LO_t> mgr;
571 solverTypesSL.sublist(FixedPoint_name).setParameters(
572 *mgr.getValidParameters()
573 );
574 }
575 {
576 Thyra::BelosTpetraGmres<Scalar,MV_t,LO_t> mgr;
577 solverTypesSL.sublist(TpetraGmres_name).setParameters(
578 *mgr.getValidParameters()
579 );
580 }
581 {
582 Thyra::BelosTpetraGmresPipeline<Scalar,MV_t,LO_t> mgr;
583 solverTypesSL.sublist(TpetraGmresPipeline_name).setParameters(
584 *mgr.getValidParameters()
585 );
586 }
587 {
588 Thyra::BelosTpetraGmresSingleReduce<Scalar,MV_t,LO_t> mgr;
589 solverTypesSL.sublist(TpetraGmresSingleReduce_name).setParameters(
590 *mgr.getValidParameters()
591 );
592 }
593 {
594 Thyra::BelosTpetraGmresSstep<Scalar,MV_t,LO_t> mgr;
595 solverTypesSL.sublist(TpetraGmresSstep_name).setParameters(
596 *mgr.getValidParameters()
597 );
598 }
599 }
600 return validParamList;
601}
602
603
604template<class Scalar>
605void BelosLinearOpWithSolveFactory<Scalar>::updateThisValidParamList()
606{
607 thisValidParamList_ = Teuchos::rcp(
608 new Teuchos::ParameterList(*generateAndGetValidParameters())
609 );
610 Teuchos::setupVerboseObjectSublist(&*thisValidParamList_);
611}
612
613
614template<class Scalar>
615void BelosLinearOpWithSolveFactory<Scalar>::initializeOpImpl(
616 const RCP<const LinearOpSourceBase<Scalar> > &fwdOpSrc,
617 const RCP<const LinearOpSourceBase<Scalar> > &approxFwdOpSrc,
618 const RCP<const PreconditionerBase<Scalar> > &prec_in,
619 const bool reusePrec,
620 LinearOpWithSolveBase<Scalar> *Op,
621 const ESupportSolveUse supportSolveUse
622 ) const
623{
624
625 using Teuchos::rcp;
626 using Teuchos::set_extra_data;
627 typedef Teuchos::ScalarTraits<Scalar> ST;
628 typedef MultiVectorBase<Scalar> MV_t;
629 typedef LinearOpBase<Scalar> LO_t;
630
631 const RCP<Teuchos::FancyOStream> out = this->getOStream();
632 const Teuchos::EVerbosityLevel verbLevel = this->getVerbLevel();
633 Teuchos::OSTab tab(out);
634 if(out.get() && static_cast<int>(verbLevel) > static_cast<int>(Teuchos::VERB_LOW))
635 *out << "\nEntering Thyra::BelosLinearOpWithSolveFactory<"<<ST::name()<<">::initializeOpImpl(...) ...\n";
636
637 // These lines are changing the verbosity of the preconditioner, which has its own verbose object list,
638 // so I am commenting these out, as it is not the job of the linear solver to dictate preconditioner verbosity.
639 //typedef Teuchos::VerboseObjectTempState<PreconditionerFactoryBase<Scalar> > VOTSPF;
640 //VOTSPF precFactoryOutputTempState(precFactory_,out,verbLevel);
641
642 TEUCHOS_TEST_FOR_EXCEPT(Op==NULL);
643 TEUCHOS_TEST_FOR_EXCEPT(fwdOpSrc.get()==NULL);
644 TEUCHOS_TEST_FOR_EXCEPT(fwdOpSrc->getOp().get()==NULL);
645 RCP<const LinearOpBase<Scalar> >
646 fwdOp = fwdOpSrc->getOp(),
647 approxFwdOp = ( approxFwdOpSrc.get() ? approxFwdOpSrc->getOp() : Teuchos::null );
648
649 //
650 // Get the BelosLinearOpWithSolve interface
651 //
652
654 *belosOp = &Teuchos::dyn_cast<BelosLinearOpWithSolve<Scalar> >(*Op);
655
656 //
657 // Get/Create the preconditioner
658 //
659
660 RCP<PreconditionerBase<Scalar> > myPrec = Teuchos::null;
661 RCP<const PreconditionerBase<Scalar> > prec = Teuchos::null;
662 if(prec_in.get()) {
663 // Use an externally defined preconditioner
664 prec = prec_in;
665 }
666 else {
667 // Try and generate a preconditioner on our own
668 if(precFactory_.get()) {
669 myPrec =
670 ( !belosOp->isExternalPrec()
671 ? Teuchos::rcp_const_cast<PreconditionerBase<Scalar> >(belosOp->extract_prec())
672 : Teuchos::null
673 );
674 bool hasExistingPrec = false;
675 if(myPrec.get()) {
676 hasExistingPrec = true;
677 // ToDo: Get the forward operator and validate that it is the same
678 // operator that is used here!
679 }
680 else {
681 hasExistingPrec = false;
682 myPrec = precFactory_->createPrec();
683 }
684 if( hasExistingPrec && reusePrec ) {
685 // Just reuse the existing preconditioner again!
686 }
687 else {
688 // Update the preconditioner
689 if(approxFwdOp.get())
690 precFactory_->initializePrec(approxFwdOpSrc,&*myPrec);
691 else
692 precFactory_->initializePrec(fwdOpSrc,&*myPrec);
693 }
694 prec = myPrec;
695 }
696 }
697
698 //
699 // Uninitialize the current solver object
700 //
701
702 bool oldIsExternalPrec = false;
703 RCP<Belos::LinearProblem<Scalar,MV_t,LO_t> > oldLP = Teuchos::null;
704 RCP<Belos::SolverManager<Scalar,MV_t,LO_t> > oldIterSolver = Teuchos::null;
705 RCP<const LinearOpSourceBase<Scalar> > oldFwdOpSrc = Teuchos::null;
706 RCP<const LinearOpSourceBase<Scalar> > oldApproxFwdOpSrc = Teuchos::null;
707 ESupportSolveUse oldSupportSolveUse = SUPPORT_SOLVE_UNSPECIFIED;
708
709 belosOp->uninitialize( &oldLP, NULL, &oldIterSolver, &oldFwdOpSrc,
710 NULL, &oldIsExternalPrec, &oldApproxFwdOpSrc, &oldSupportSolveUse );
711
712 //
713 // Create the Belos linear problem
714 // NOTE: If one exists already, reuse it.
715 //
716
717 typedef Belos::LinearProblem<Scalar,MV_t,LO_t> LP_t;
718 RCP<LP_t> lp;
719 if (oldLP != Teuchos::null) {
720 lp = oldLP;
721 }
722 else {
723 lp = rcp(new LP_t());
724 }
725
726 //
727 // Set the operator
728 //
729
730 lp->setOperator(fwdOp);
731
732 //
733 // Set the preconditioner
734 //
735
736 if(prec.get()) {
737 RCP<const LinearOpBase<Scalar> > unspecified = prec->getUnspecifiedPrecOp();
738 RCP<const LinearOpBase<Scalar> > left = prec->getLeftPrecOp();
739 RCP<const LinearOpBase<Scalar> > right = prec->getRightPrecOp();
740 TEUCHOS_TEST_FOR_EXCEPTION(
741 !( left.get() || right.get() || unspecified.get() ), std::logic_error
742 ,"Error, at least one preconditoner linear operator objects must be set!"
743 );
744 if (nonnull(unspecified)) {
745 if (paramList_->get<bool>(LeftPreconditionerIfUnspecified_name, false))
746 lp->setLeftPrec(unspecified);
747 else
748 lp->setRightPrec(unspecified);
749 }
750 else if (nonnull(left)) {
751 lp->setLeftPrec(left);
752 }
753 else if (nonnull(right)) {
754 lp->setRightPrec(right);
755 }
756 else {
757 // Set a left, right or split preconditioner
758 TEUCHOS_TEST_FOR_EXCEPTION(
759 nonnull(left) && nonnull(right),std::logic_error
760 ,"Error, we can not currently handle split preconditioners!"
761 );
762 }
763 }
764 if(myPrec.get()) {
765 set_extra_data<RCP<PreconditionerBase<Scalar> > >(myPrec,"Belos::InternalPrec",
766 Teuchos::inOutArg(lp), Teuchos::POST_DESTROY, false);
767 }
768 else if(prec.get()) {
769 set_extra_data<RCP<const PreconditionerBase<Scalar> > >(prec,"Belos::ExternalPrec",
770 Teuchos::inOutArg(lp), Teuchos::POST_DESTROY, false);
771 }
772
773 //
774 // Generate the parameter list.
775 //
776
777 typedef Belos::SolverManager<Scalar,MV_t,LO_t> IterativeSolver_t;
778 RCP<IterativeSolver_t> iterativeSolver = Teuchos::null;
779 RCP<Teuchos::ParameterList> solverPL = Teuchos::rcp( new Teuchos::ParameterList() );
780
781 switch(solverType_) {
782 case SOLVER_TYPE_BLOCK_GMRES:
783 {
784 // Set the PL
785 if(paramList_.get()) {
786 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
787 Teuchos::ParameterList &gmresPL = solverTypesPL.sublist(BlockGMRES_name);
788 solverPL = Teuchos::rcp( &gmresPL, false );
789 }
790 // Create the solver
791 if (oldIterSolver != Teuchos::null) {
792 iterativeSolver = oldIterSolver;
793 iterativeSolver->setProblem( lp );
794 iterativeSolver->setParameters( solverPL );
795 }
796 else {
797 iterativeSolver = rcp(new Belos::BlockGmresSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
798 }
799 break;
800 }
801 case SOLVER_TYPE_PSEUDO_BLOCK_GMRES:
802 {
803 // Set the PL
804 if(paramList_.get()) {
805 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
806 Teuchos::ParameterList &pbgmresPL = solverTypesPL.sublist(PseudoBlockGMRES_name);
807 solverPL = Teuchos::rcp( &pbgmresPL, false );
808 }
809 //
810 // Create the solver
811 //
812 if (oldIterSolver != Teuchos::null) {
813 iterativeSolver = oldIterSolver;
814 iterativeSolver->setProblem( lp );
815 iterativeSolver->setParameters( solverPL );
816 }
817 else {
818 iterativeSolver = rcp(new Belos::PseudoBlockGmresSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
819 }
820 break;
821 }
822 case SOLVER_TYPE_BLOCK_CG:
823 {
824 // Set the PL
825 if(paramList_.get()) {
826 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
827 Teuchos::ParameterList &cgPL = solverTypesPL.sublist(BlockCG_name);
828 solverPL = Teuchos::rcp( &cgPL, false );
829 }
830 // Create the solver
831 if (oldIterSolver != Teuchos::null) {
832 iterativeSolver = oldIterSolver;
833 iterativeSolver->setProblem( lp );
834 iterativeSolver->setParameters( solverPL );
835 }
836 else {
837 iterativeSolver = rcp(new Belos::BlockCGSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
838 }
839 break;
840 }
841 case SOLVER_TYPE_PSEUDO_BLOCK_CG:
842 {
843 // Set the PL
844 if(paramList_.get()) {
845 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
846 Teuchos::ParameterList &pbcgPL = solverTypesPL.sublist(PseudoBlockCG_name);
847 solverPL = Teuchos::rcp( &pbcgPL, false );
848 }
849 //
850 // Create the solver
851 //
852 if (oldIterSolver != Teuchos::null) {
853 iterativeSolver = oldIterSolver;
854 iterativeSolver->setProblem( lp );
855 iterativeSolver->setParameters( solverPL );
856 }
857 else {
858 iterativeSolver = rcp(new Belos::PseudoBlockCGSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
859 }
860 break;
861 }
862 case SOLVER_TYPE_PSEUDO_BLOCK_STOCHASTIC_CG:
863 {
864 // Set the PL
865 if(paramList_.get()) {
866 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
867 Teuchos::ParameterList &pbcgPL = solverTypesPL.sublist(PseudoBlockStochasticCG_name);
868 solverPL = Teuchos::rcp( &pbcgPL, false );
869 }
870 //
871 // Create the solver
872 //
873 if (oldIterSolver != Teuchos::null) {
874 iterativeSolver = oldIterSolver;
875 iterativeSolver->setProblem( lp );
876 iterativeSolver->setParameters( solverPL );
877 }
878 else {
879 iterativeSolver = rcp(new Belos::PseudoBlockStochasticCGSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
880 }
881 break;
882 }
883 case SOLVER_TYPE_GCRODR:
884 {
885 // Set the PL
886 if(paramList_.get()) {
887 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
888 Teuchos::ParameterList &gcrodrPL = solverTypesPL.sublist(GCRODR_name);
889 solverPL = Teuchos::rcp( &gcrodrPL, false );
890 }
891 // Create the solver
892 if (oldIterSolver != Teuchos::null) {
893 iterativeSolver = oldIterSolver;
894 iterativeSolver->setProblem( lp );
895 iterativeSolver->setParameters( solverPL );
896 }
897 else {
898 iterativeSolver = rcp(new Belos::GCRODRSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
899 }
900 break;
901 }
902 case SOLVER_TYPE_RCG:
903 {
904 // Set the PL
905 if(paramList_.get()) {
906 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
907 Teuchos::ParameterList &rcgPL = solverTypesPL.sublist(RCG_name);
908 solverPL = Teuchos::rcp( &rcgPL, false );
909 }
910 // Create the solver
911 if (oldIterSolver != Teuchos::null) {
912 iterativeSolver = oldIterSolver;
913 iterativeSolver->setProblem( lp );
914 iterativeSolver->setParameters( solverPL );
915 }
916 else {
917 iterativeSolver = rcp(new Belos::RCGSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
918 }
919 break;
920 }
921 case SOLVER_TYPE_MINRES:
922 {
923 // Set the PL
924 if(paramList_.get()) {
925 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
926 Teuchos::ParameterList &minresPL = solverTypesPL.sublist(MINRES_name);
927 solverPL = Teuchos::rcp( &minresPL, false );
928 }
929 // Create the solver
930 if (oldIterSolver != Teuchos::null) {
931 iterativeSolver = oldIterSolver;
932 iterativeSolver->setProblem( lp );
933 iterativeSolver->setParameters( solverPL );
934 }
935 else {
936 iterativeSolver = rcp(new Belos::MinresSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
937 }
938 break;
939 }
940 case SOLVER_TYPE_TFQMR:
941 {
942 // Set the PL
943 if(paramList_.get()) {
944 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
945 Teuchos::ParameterList &minresPL = solverTypesPL.sublist(TFQMR_name);
946 solverPL = Teuchos::rcp( &minresPL, false );
947 }
948 // Create the solver
949 if (oldIterSolver != Teuchos::null) {
950 iterativeSolver = oldIterSolver;
951 iterativeSolver->setProblem( lp );
952 iterativeSolver->setParameters( solverPL );
953 }
954 else {
955 iterativeSolver = rcp(new Belos::TFQMRSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
956 }
957 break;
958 }
959 case SOLVER_TYPE_BICGSTAB:
960 {
961 // Set the PL
962 if(paramList_.get()) {
963 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
964 Teuchos::ParameterList &bicgstabPL = solverTypesPL.sublist(BiCGStab_name);
965 solverPL = Teuchos::rcp( &bicgstabPL, false );
966 }
967 // Create the solver
968 if (oldIterSolver != Teuchos::null) {
969 iterativeSolver = oldIterSolver;
970 iterativeSolver->setProblem( lp );
971 iterativeSolver->setParameters( solverPL );
972 }
973 else {
974 iterativeSolver = rcp(new Belos::BiCGStabSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
975 }
976 break;
977 }
978 case SOLVER_TYPE_FIXEDPOINT:
979 {
980 // Set the PL
981 if(paramList_.get()) {
982 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
983 Teuchos::ParameterList &fixedPointPL = solverTypesPL.sublist(FixedPoint_name);
984 solverPL = Teuchos::rcp( &fixedPointPL, false );
985 }
986 // Create the solver
987 if (oldIterSolver != Teuchos::null) {
988 iterativeSolver = oldIterSolver;
989 iterativeSolver->setProblem( lp );
990 iterativeSolver->setParameters( solverPL );
991 }
992 else {
993 iterativeSolver = rcp(new Belos::FixedPointSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
994 }
995 break;
996 }
997 case SOLVER_TYPE_TPETRA_GMRES:
998 {
999 // Get the PL
1000 if(paramList_.get()) {
1001 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
1002 Teuchos::ParameterList &tpetraGmresPL = solverTypesPL.sublist(TpetraGmres_name);
1003 solverPL = Teuchos::rcp( &tpetraGmresPL, false );
1004 }
1005 // Create the solver, always
1006 iterativeSolver = rcp(new Thyra::BelosTpetraGmres<Scalar,MV_t,LO_t>());
1007 // Set the Problem & PL
1008 iterativeSolver->setProblem( lp );
1009 iterativeSolver->setParameters( solverPL );
1010 break;
1011 }
1012 case SOLVER_TYPE_TPETRA_GMRES_PIPELINE:
1013 {
1014 // Get the PL
1015 if(paramList_.get()) {
1016 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
1017 Teuchos::ParameterList &tpetraGmresPipelinePL = solverTypesPL.sublist(TpetraGmresPipeline_name);
1018 solverPL = Teuchos::rcp( &tpetraGmresPipelinePL, false );
1019 }
1020 // Create the solver, always
1021 iterativeSolver = rcp(new Thyra::BelosTpetraGmresPipeline<Scalar,MV_t,LO_t>());
1022 // Set the Problem & PL
1023 iterativeSolver->setProblem( lp );
1024 iterativeSolver->setParameters( solverPL );
1025 break;
1026 }
1027 case SOLVER_TYPE_TPETRA_GMRES_SINGLE_REDUCE:
1028 {
1029 // Get the PL
1030 if(paramList_.get()) {
1031 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
1032 Teuchos::ParameterList &tpetraGmresSingleReducePL = solverTypesPL.sublist(TpetraGmresSingleReduce_name);
1033 solverPL = Teuchos::rcp( &tpetraGmresSingleReducePL, false );
1034 }
1035 // Create the solver, always
1036 iterativeSolver = rcp(new Thyra::BelosTpetraGmresSingleReduce<Scalar,MV_t,LO_t>());
1037 // Set the Problem & PL
1038 iterativeSolver->setProblem( lp );
1039 iterativeSolver->setParameters( solverPL );
1040 break;
1041 }
1042 case SOLVER_TYPE_TPETRA_GMRES_SSTEP:
1043 {
1044 // Get the PL
1045 if(paramList_.get()) {
1046 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
1047 Teuchos::ParameterList &tpetraGmresSstepPL = solverTypesPL.sublist(TpetraGmresSstep_name);
1048 solverPL = Teuchos::rcp( &tpetraGmresSstepPL, false );
1049 }
1050 // Create the solver, always
1051 iterativeSolver = rcp(new Thyra::BelosTpetraGmresSstep<Scalar,MV_t,LO_t>());
1052 // Set the Problem & PL
1053 iterativeSolver->setProblem( lp );
1054 iterativeSolver->setParameters( solverPL );
1055 break;
1056 }
1057
1058 default:
1059 {
1060 TEUCHOS_TEST_FOR_EXCEPT(true);
1061 }
1062 }
1063
1064 //
1065 // Initialize the LOWS object
1066 //
1067
1068 belosOp->initialize(
1069 lp, solverPL, iterativeSolver,
1070 fwdOpSrc, prec, myPrec.get()==NULL, approxFwdOpSrc,
1071 supportSolveUse, convergenceTestFrequency_
1072 );
1073 belosOp->setOStream(out);
1074 belosOp->setVerbLevel(verbLevel);
1075#ifdef TEUCHOS_DEBUG
1076 if(paramList_.get()) {
1077 // Make sure we read the list correctly
1078 paramList_->validateParameters(*this->getValidParameters(),1); // Validate 0th and 1st level deep
1079 }
1080#endif
1081 if(out.get() && static_cast<int>(verbLevel) > static_cast<int>(Teuchos::VERB_LOW))
1082 *out << "\nLeaving Thyra::BelosLinearOpWithSolveFactory<"<<ST::name()<<">::initializeOpImpl(...) ...\n";
1083
1084}
1085
1086
1087} // namespace Thyra
1088
1089
1090#endif // THYRA_BELOS_LINEAR_OP_WITH_SOLVE_FACTORY_HPP
Thyra specializations of MultiVecTraits and OperatorTraits.
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > &params) override
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > &params) override
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > &params) override
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > &params) override
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > &params) override
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > &params) override
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > &params) override
void unsetPreconditionerFactory(Teuchos::RCP< PreconditionerFactoryBase< Scalar > > *precFactory, std::string *precFactoryName)
bool supportsPreconditionerInputType(const EPreconditionerInputType precOpType) const
bool isCompatible(const LinearOpSourceBase< Scalar > &fwdOpSrc) const
void setParameterList(Teuchos::RCP< Teuchos::ParameterList > const &paramList)
void initializePreconditionedOp(const Teuchos::RCP< const LinearOpSourceBase< Scalar > > &fwdOpSrc, const Teuchos::RCP< const PreconditionerBase< Scalar > > &prec, LinearOpWithSolveBase< Scalar > *Op, const ESupportSolveUse supportSolveUse) const
void setPreconditionerFactory(const Teuchos::RCP< PreconditionerFactoryBase< Scalar > > &precFactory, const std::string &precFactoryName)
void initializeAndReuseOp(const Teuchos::RCP< const LinearOpSourceBase< Scalar > > &fwdOpSrc, LinearOpWithSolveBase< Scalar > *Op) const
Teuchos::RCP< const Teuchos::ParameterList > getParameterList() const
void initializeOp(const Teuchos::RCP< const LinearOpSourceBase< Scalar > > &fwdOpSrc, LinearOpWithSolveBase< Scalar > *Op, const ESupportSolveUse supportSolveUse) const
Teuchos::RCP< Teuchos::ParameterList > getNonconstParameterList()
void uninitializeOp(LinearOpWithSolveBase< Scalar > *Op, Teuchos::RCP< const LinearOpSourceBase< Scalar > > *fwdOpSrc, Teuchos::RCP< const PreconditionerBase< Scalar > > *prec, Teuchos::RCP< const LinearOpSourceBase< Scalar > > *approxFwdOpSrc, ESupportSolveUse *supportSolveUse) const
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
BelosLinearOpWithSolveFactory()
Construct without preconditioner factory.
void initializeApproxPreconditionedOp(const Teuchos::RCP< const LinearOpSourceBase< Scalar > > &fwdOpSrc, const Teuchos::RCP< const LinearOpSourceBase< Scalar > > &approxFwdOpSrc, LinearOpWithSolveBase< Scalar > *Op, const ESupportSolveUse supportSolveUse) const
Teuchos::RCP< LinearOpWithSolveBase< Scalar > > createOp() const
Teuchos::RCP< Teuchos::ParameterList > unsetParameterList()
Teuchos::RCP< PreconditionerFactoryBase< Scalar > > getPreconditionerFactory() const
Concrete LinearOpWithSolveBase subclass in terms of Belos.
void uninitialize(RCP< Belos::LinearProblem< Scalar, MV_t, LO_t > > *lp=NULL, RCP< Teuchos::ParameterList > *solverPL=NULL, RCP< Belos::SolverManager< Scalar, MV_t, LO_t > > *iterativeSolver=NULL, RCP< const LinearOpSourceBase< Scalar > > *fwdOpSrc=NULL, RCP< const PreconditionerBase< Scalar > > *prec=NULL, bool *isExternalPrec=NULL, RCP< const LinearOpSourceBase< Scalar > > *approxFwdOpSrc=NULL, ESupportSolveUse *supportSolveUse=NULL)
Uninitializes and returns stored quantities.
RCP< const LinearOpSourceBase< Scalar > > extract_approxFwdOpSrc()
RCP< const PreconditionerBase< Scalar > > extract_prec()
RCP< const LinearOpSourceBase< Scalar > > extract_fwdOpSrc()

Generated for Stratimikos by doxygen 1.17.0