Stratimikos Version of the Day
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
65#include "BelosThyraAdapter.hpp"
66#include "Teuchos_VerboseObjectParameterListHelpers.hpp"
67#include "Teuchos_StandardParameterEntryValidators.hpp"
68#include "Teuchos_ParameterList.hpp"
69#include "Teuchos_dyn_cast.hpp"
70#include "Teuchos_ValidatorXMLConverterDB.hpp"
71#include "Teuchos_StandardValidatorXMLConverters.hpp"
72
73
74namespace Thyra {
75
76
77// Parameter names for Parameter List
78
79template<class Scalar>
80const std::string BelosLinearOpWithSolveFactory<Scalar>::SolverType_name = "Solver Type";
81template<class Scalar>
82const std::string BelosLinearOpWithSolveFactory<Scalar>::SolverType_default = "Pseudo Block GMRES";
83template<class Scalar>
84const std::string BelosLinearOpWithSolveFactory<Scalar>::SolverTypes_name = "Solver Types";
85template<class Scalar>
86const std::string BelosLinearOpWithSolveFactory<Scalar>::BlockGMRES_name = "Block GMRES";
87template<class Scalar>
88const std::string BelosLinearOpWithSolveFactory<Scalar>::PseudoBlockGMRES_name = "Pseudo Block GMRES";
89template<class Scalar>
90const std::string BelosLinearOpWithSolveFactory<Scalar>::BlockCG_name = "Block CG";
91template<class Scalar>
92const std::string BelosLinearOpWithSolveFactory<Scalar>::PseudoBlockCG_name = "Pseudo Block CG";
93template<class Scalar>
94const std::string BelosLinearOpWithSolveFactory<Scalar>::PseudoBlockStochasticCG_name = "Pseudo Block Stochastic CG";
95template<class Scalar>
97template<class Scalar>
99template<class Scalar>
100const std::string BelosLinearOpWithSolveFactory<Scalar>::MINRES_name = "MINRES";
101template<class Scalar>
102const std::string BelosLinearOpWithSolveFactory<Scalar>::TFQMR_name = "TFQMR";
103template<class Scalar>
104const std::string BelosLinearOpWithSolveFactory<Scalar>::BiCGStab_name = "BiCGStab";
105template<class Scalar>
106const std::string BelosLinearOpWithSolveFactory<Scalar>::FixedPoint_name = "Fixed Point";
107template<class Scalar>
108const std::string BelosLinearOpWithSolveFactory<Scalar>::ConvergenceTestFrequency_name = "Convergence Test Frequency";
109
110namespace {
111const std::string LeftPreconditionerIfUnspecified_name = "Left Preconditioner If Unspecified";
112}
113
114// Constructors/initializers/accessors
115
116
117template<class Scalar>
119 :solverType_(SOLVER_TYPE_PSEUDO_BLOCK_GMRES),
120 convergenceTestFrequency_(1)
121{
122 updateThisValidParamList();
123}
124
125
126template<class Scalar>
128 const RCP<PreconditionerFactoryBase<Scalar> > &precFactory
129 )
130 :solverType_(SOLVER_TYPE_PSEUDO_BLOCK_GMRES)
131{
132 this->setPreconditionerFactory(precFactory, "");
133}
134
135
136// Overridden from LinearOpWithSolveFactoryBase
137
138
139template<class Scalar>
141{
142 return true;
143}
144
145
146template<class Scalar>
148 const RCP<PreconditionerFactoryBase<Scalar> > &precFactory,
149 const std::string &precFactoryName
150 )
151{
152 TEUCHOS_TEST_FOR_EXCEPT(!precFactory.get());
153 RCP<const Teuchos::ParameterList>
154 precFactoryValidPL = precFactory->getValidParameters();
155 const std::string _precFactoryName =
156 ( precFactoryName != ""
157 ? precFactoryName
158 : ( precFactoryValidPL.get() ? precFactoryValidPL->name() : "GENERIC PRECONDITIONER FACTORY" )
159 );
160 precFactory_ = precFactory;
161 precFactoryName_ = _precFactoryName;
162 updateThisValidParamList();
163}
164
165
166template<class Scalar>
167RCP<PreconditionerFactoryBase<Scalar> >
169{
170 return precFactory_;
171}
172
173
174template<class Scalar>
176 RCP<PreconditionerFactoryBase<Scalar> > *precFactory,
177 std::string *precFactoryName
178 )
179{
180 if(precFactory) *precFactory = precFactory_;
181 if(precFactoryName) *precFactoryName = precFactoryName_;
182 precFactory_ = Teuchos::null;
183 precFactoryName_ = "";
184 updateThisValidParamList();
185}
186
187
188template<class Scalar>
190 const LinearOpSourceBase<Scalar> &fwdOpSrc
191 ) const
192{
193 if(precFactory_.get())
194 return precFactory_->isCompatible(fwdOpSrc);
195 return true; // Without a preconditioner, we are compatible with all linear operators!
196}
197
198
199template<class Scalar>
200RCP<LinearOpWithSolveBase<Scalar> >
202{
203 return Teuchos::rcp(new BelosLinearOpWithSolve<Scalar>());
204}
205
206
207template<class Scalar>
209 const RCP<const LinearOpSourceBase<Scalar> > &fwdOpSrc,
210 LinearOpWithSolveBase<Scalar> *Op,
211 const ESupportSolveUse supportSolveUse
212 ) const
213{
214 using Teuchos::null;
215 initializeOpImpl(fwdOpSrc,null,null,false,Op,supportSolveUse);
216}
217
218
219template<class Scalar>
221 const RCP<const LinearOpSourceBase<Scalar> > &fwdOpSrc,
222 LinearOpWithSolveBase<Scalar> *Op
223 ) const
224{
225 using Teuchos::null;
226 initializeOpImpl(fwdOpSrc,null,null,true,Op,SUPPORT_SOLVE_UNSPECIFIED);
227}
228
229
230template<class Scalar>
232 const EPreconditionerInputType precOpType
233 ) const
234{
235 if(precFactory_.get())
236 return true;
237 return (precOpType==PRECONDITIONER_INPUT_TYPE_AS_OPERATOR);
238}
239
240
241template<class Scalar>
243 const RCP<const LinearOpSourceBase<Scalar> > &fwdOpSrc,
244 const RCP<const PreconditionerBase<Scalar> > &prec,
245 LinearOpWithSolveBase<Scalar> *Op,
246 const ESupportSolveUse supportSolveUse
247 ) const
248{
249 using Teuchos::null;
250 initializeOpImpl(fwdOpSrc,null,prec,false,Op,supportSolveUse);
251}
252
253
254template<class Scalar>
256 const RCP<const LinearOpSourceBase<Scalar> > &fwdOpSrc,
257 const RCP<const LinearOpSourceBase<Scalar> > &approxFwdOpSrc,
258 LinearOpWithSolveBase<Scalar> *Op,
259 const ESupportSolveUse supportSolveUse
260 ) const
261{
262 using Teuchos::null;
263 initializeOpImpl(fwdOpSrc,approxFwdOpSrc,null,false,Op,supportSolveUse);
264}
265
266
267template<class Scalar>
269 LinearOpWithSolveBase<Scalar> *Op,
270 RCP<const LinearOpSourceBase<Scalar> > *fwdOpSrc,
271 RCP<const PreconditionerBase<Scalar> > *prec,
272 RCP<const LinearOpSourceBase<Scalar> > *approxFwdOpSrc,
273 ESupportSolveUse *supportSolveUse
274 ) const
275{
276#ifdef TEUCHOS_DEBUG
277 TEUCHOS_TEST_FOR_EXCEPT(Op==NULL);
278#endif
280 &belosOp = Teuchos::dyn_cast<BelosLinearOpWithSolve<Scalar> >(*Op);
281 RCP<const LinearOpSourceBase<Scalar> >
282 _fwdOpSrc = belosOp.extract_fwdOpSrc();
283 RCP<const PreconditionerBase<Scalar> >
284 _prec = ( belosOp.isExternalPrec() ? belosOp.extract_prec() : Teuchos::null );
285 // Note: above we only extract the preconditioner if it was passed in
286 // externally. Otherwise, we need to hold on to it so that we can reuse it
287 // in the next initialization.
288 RCP<const LinearOpSourceBase<Scalar> >
289 _approxFwdOpSrc = belosOp.extract_approxFwdOpSrc();
290 ESupportSolveUse
291 _supportSolveUse = belosOp.supportSolveUse();
292 if(fwdOpSrc) *fwdOpSrc = _fwdOpSrc;
293 if(prec) *prec = _prec;
294 if(approxFwdOpSrc) *approxFwdOpSrc = _approxFwdOpSrc;
295 if(supportSolveUse) *supportSolveUse = _supportSolveUse;
296}
297
298
299// Overridden from ParameterListAcceptor
300
301
302template<class Scalar>
304 RCP<Teuchos::ParameterList> const& paramList
305 )
306{
307 TEUCHOS_TEST_FOR_EXCEPT(paramList.get()==NULL);
308 paramList->validateParametersAndSetDefaults(*this->getValidParameters(), 1);
309 paramList_ = paramList;
310 solverType_ =
311 Teuchos::getIntegralValue<EBelosSolverType>(*paramList_, SolverType_name);
312 convergenceTestFrequency_ =
313 Teuchos::getParameter<int>(*paramList_, ConvergenceTestFrequency_name);
314 Teuchos::readVerboseObjectSublist(&*paramList_,this);
315}
316
317
318template<class Scalar>
319RCP<Teuchos::ParameterList>
321{
322 return paramList_;
323}
324
325
326template<class Scalar>
327RCP<Teuchos::ParameterList>
329{
330 RCP<Teuchos::ParameterList> _paramList = paramList_;
331 paramList_ = Teuchos::null;
332 return _paramList;
333}
334
335
336template<class Scalar>
337RCP<const Teuchos::ParameterList>
339{
340 return paramList_;
341}
342
343
344template<class Scalar>
345RCP<const Teuchos::ParameterList>
347{
348 return thisValidParamList_;
349}
350
351
352// Public functions overridden from Teuchos::Describable
353
354
355template<class Scalar>
357{
358 std::ostringstream oss;
359 oss << "Thyra::BelosLinearOpWithSolveFactory";
360 //oss << "{";
361 // ToDo: Fill this in some!
362 //oss << "}";
363 return oss.str();
364}
365
366
367// private
368
369
370template<class Scalar>
371RCP<const Teuchos::ParameterList>
373{
374 using Teuchos::as;
375 using Teuchos::tuple;
376 using Teuchos::setStringToIntegralParameter;
377Teuchos::ValidatorXMLConverterDB::addConverter(
378 Teuchos::DummyObjectGetter<
379 Teuchos::StringToIntegralParameterEntryValidator<EBelosSolverType>
380 >::getDummyObject(),
381 Teuchos::DummyObjectGetter<Teuchos::StringToIntegralValidatorXMLConverter<
382 EBelosSolverType> >::getDummyObject());
383
384 typedef MultiVectorBase<Scalar> MV_t;
385 typedef LinearOpBase<Scalar> LO_t;
386 static RCP<Teuchos::ParameterList> validParamList;
387 if(validParamList.get()==NULL) {
388 validParamList = Teuchos::rcp(new Teuchos::ParameterList("BelosLinearOpWithSolveFactory"));
389 setStringToIntegralParameter<EBelosSolverType>(
390 SolverType_name, SolverType_default,
391 "Type of linear solver algorithm to use.",
392 tuple<std::string>(
393 "Block GMRES",
394 "Pseudo Block GMRES",
395 "Block CG",
396 "Pseudo Block CG",
397 "Pseudo Block Stochastic CG",
398 "GCRODR",
399 "RCG",
400 "MINRES",
401 "TFQMR",
402 "BiCGStab",
403 "Fixed Point"
404 ),
405 tuple<std::string>(
406 "Block GMRES solver for nonsymmetric linear systems. It can also solve "
407 "single right-hand side systems, and can also perform Flexible GMRES "
408 "(where the preconditioner may change at every iteration, for example "
409 "for inner-outer iterations), by setting options in the \"Block GMRES\" "
410 "sublist.",
411
412 "GMRES solver for nonsymmetric linear systems, that performs single "
413 "right-hand side solves on multiple right-hand sides at once. It "
414 "exploits operator multivector multiplication in order to amortize "
415 "global communication costs. Individual linear systems are deflated "
416 "out as they are solved.",
417
418 "Block CG solver for symmetric (Hermitian in complex arithmetic) "
419 "positive definite linear systems. It can also solve single "
420 "right-hand-side systems.",
421
422 "CG solver that performs single right-hand side CG on multiple right-hand "
423 "sides at once. It exploits operator multivector multiplication in order "
424 "to amortize global communication costs. Individual linear systems are "
425 "deflated out as they are solved.",
426
427 "stochastic CG solver that performs single right-hand side CG on multiple right-hand "
428 "sides at once. It exploits operator multivector multiplication in order "
429 "to amortize global communication costs. Individual linear systems are "
430 "deflated out as they are solved. [EXPERIMENTAL]",
431
432 "Variant of GMRES that performs subspace recycling to accelerate "
433 "convergence for sequences of solves with related linear systems. "
434 "Individual linear systems are deflated out as they are solved. "
435 "The current implementation only supports real-valued Scalar types.",
436
437 "CG solver for symmetric (Hermitian in complex arithmetic) positive "
438 "definite linear systems, that performs subspace recycling to "
439 "accelerate convergence for sequences of related linear systems.",
440
441 "MINRES solver for symmetric indefinite linear systems. It performs "
442 "single-right-hand-side solves on multiple right-hand sides sequentially.",
443
444 "TFQMR (Transpose-Free QMR) solver for nonsymmetric linear systems.",
445
446 "BiCGStab solver for nonsymmetric linear systems.",
447
448 "Fixed point iteration"
449 ),
450 tuple<EBelosSolverType>(
451 SOLVER_TYPE_BLOCK_GMRES,
452 SOLVER_TYPE_PSEUDO_BLOCK_GMRES,
453 SOLVER_TYPE_BLOCK_CG,
454 SOLVER_TYPE_PSEUDO_BLOCK_CG,
455 SOLVER_TYPE_PSEUDO_BLOCK_STOCHASTIC_CG,
456 SOLVER_TYPE_GCRODR,
457 SOLVER_TYPE_RCG,
458 SOLVER_TYPE_MINRES,
459 SOLVER_TYPE_TFQMR,
460 SOLVER_TYPE_BICGSTAB,
461 SOLVER_TYPE_FIXEDPOINT
462 ),
463 &*validParamList
464 );
465 validParamList->set(ConvergenceTestFrequency_name, as<int>(1),
466 "Number of linear solver iterations to skip between applying"
467 " user-defined convergence test.");
468 validParamList->set(
469 LeftPreconditionerIfUnspecified_name, false,
470 "If the preconditioner does not specify if it is left or right, and this\n"
471 "option is set to true, put the preconditioner on the left side.\n"
472 "Historically, preconditioning is on the right. Some solvers may not\n"
473 "support left preconditioning.");
474 Teuchos::ParameterList
475 &solverTypesSL = validParamList->sublist(SolverTypes_name);
476 {
478 solverTypesSL.sublist(BlockGMRES_name).setParameters(
479 *mgr.getValidParameters()
480 );
481 }
482 {
484 solverTypesSL.sublist(PseudoBlockGMRES_name).setParameters(
485 *mgr.getValidParameters()
486 );
487 }
488 {
490 solverTypesSL.sublist(BlockCG_name).setParameters(
491 *mgr.getValidParameters()
492 );
493 }
494 {
496 solverTypesSL.sublist(PseudoBlockCG_name).setParameters(
497 *mgr.getValidParameters()
498 );
499 }
500 {
502 solverTypesSL.sublist(PseudoBlockStochasticCG_name).setParameters(
503 *mgr.getValidParameters()
504 );
505 }
506 {
508 solverTypesSL.sublist(GCRODR_name).setParameters(
509 *mgr.getValidParameters()
510 );
511 }
512 {
514 solverTypesSL.sublist(RCG_name).setParameters(
515 *mgr.getValidParameters()
516 );
517 }
518 {
520 solverTypesSL.sublist(MINRES_name).setParameters(
521 *mgr.getValidParameters()
522 );
523 }
524 {
526 solverTypesSL.sublist(TFQMR_name).setParameters(
527 *mgr.getValidParameters()
528 );
529 }
530 {
532 solverTypesSL.sublist(BiCGStab_name).setParameters(
533 *mgr.getValidParameters()
534 );
535 }
536 {
538 solverTypesSL.sublist(FixedPoint_name).setParameters(
539 *mgr.getValidParameters()
540 );
541 }
542 }
543 return validParamList;
544}
545
546
547template<class Scalar>
548void BelosLinearOpWithSolveFactory<Scalar>::updateThisValidParamList()
549{
550 thisValidParamList_ = Teuchos::rcp(
551 new Teuchos::ParameterList(*generateAndGetValidParameters())
552 );
553 Teuchos::setupVerboseObjectSublist(&*thisValidParamList_);
554}
555
556
557template<class Scalar>
558void BelosLinearOpWithSolveFactory<Scalar>::initializeOpImpl(
559 const RCP<const LinearOpSourceBase<Scalar> > &fwdOpSrc,
560 const RCP<const LinearOpSourceBase<Scalar> > &approxFwdOpSrc,
561 const RCP<const PreconditionerBase<Scalar> > &prec_in,
562 const bool reusePrec,
563 LinearOpWithSolveBase<Scalar> *Op,
564 const ESupportSolveUse supportSolveUse
565 ) const
566{
567
568 using Teuchos::rcp;
569 using Teuchos::set_extra_data;
570 typedef Teuchos::ScalarTraits<Scalar> ST;
571 typedef MultiVectorBase<Scalar> MV_t;
572 typedef LinearOpBase<Scalar> LO_t;
573
574 const RCP<Teuchos::FancyOStream> out = this->getOStream();
575 const Teuchos::EVerbosityLevel verbLevel = this->getVerbLevel();
576 Teuchos::OSTab tab(out);
577 if(out.get() && static_cast<int>(verbLevel) > static_cast<int>(Teuchos::VERB_LOW))
578 *out << "\nEntering Thyra::BelosLinearOpWithSolveFactory<"<<ST::name()<<">::initializeOpImpl(...) ...\n";
579
580 // These lines are changing the verbosity of the preconditioner, which has its own verbose object list,
581 // so I am commenting these out, as it is not the job of the linear solver to dictate preconditioner verbosity.
582 //typedef Teuchos::VerboseObjectTempState<PreconditionerFactoryBase<Scalar> > VOTSPF;
583 //VOTSPF precFactoryOutputTempState(precFactory_,out,verbLevel);
584
585 TEUCHOS_TEST_FOR_EXCEPT(Op==NULL);
586 TEUCHOS_TEST_FOR_EXCEPT(fwdOpSrc.get()==NULL);
587 TEUCHOS_TEST_FOR_EXCEPT(fwdOpSrc->getOp().get()==NULL);
588 RCP<const LinearOpBase<Scalar> >
589 fwdOp = fwdOpSrc->getOp(),
590 approxFwdOp = ( approxFwdOpSrc.get() ? approxFwdOpSrc->getOp() : Teuchos::null );
591
592 //
593 // Get the BelosLinearOpWithSolve interface
594 //
595
596 BelosLinearOpWithSolve<Scalar>
597 *belosOp = &Teuchos::dyn_cast<BelosLinearOpWithSolve<Scalar> >(*Op);
598
599 //
600 // Get/Create the preconditioner
601 //
602
603 RCP<PreconditionerBase<Scalar> > myPrec = Teuchos::null;
604 RCP<const PreconditionerBase<Scalar> > prec = Teuchos::null;
605 if(prec_in.get()) {
606 // Use an externally defined preconditioner
607 prec = prec_in;
608 }
609 else {
610 // Try and generate a preconditioner on our own
611 if(precFactory_.get()) {
612 myPrec =
613 ( !belosOp->isExternalPrec()
614 ? Teuchos::rcp_const_cast<PreconditionerBase<Scalar> >(belosOp->extract_prec())
615 : Teuchos::null
616 );
617 bool hasExistingPrec = false;
618 if(myPrec.get()) {
619 hasExistingPrec = true;
620 // ToDo: Get the forward operator and validate that it is the same
621 // operator that is used here!
622 }
623 else {
624 hasExistingPrec = false;
625 myPrec = precFactory_->createPrec();
626 }
627 if( hasExistingPrec && reusePrec ) {
628 // Just reuse the existing preconditioner again!
629 }
630 else {
631 // Update the preconditioner
632 if(approxFwdOp.get())
633 precFactory_->initializePrec(approxFwdOpSrc,&*myPrec);
634 else
635 precFactory_->initializePrec(fwdOpSrc,&*myPrec);
636 }
637 prec = myPrec;
638 }
639 }
640
641 //
642 // Uninitialize the current solver object
643 //
644
645 bool oldIsExternalPrec = false;
646 RCP<Belos::LinearProblem<Scalar,MV_t,LO_t> > oldLP = Teuchos::null;
647 RCP<Belos::SolverManager<Scalar,MV_t,LO_t> > oldIterSolver = Teuchos::null;
648 RCP<const LinearOpSourceBase<Scalar> > oldFwdOpSrc = Teuchos::null;
649 RCP<const LinearOpSourceBase<Scalar> > oldApproxFwdOpSrc = Teuchos::null;
650 ESupportSolveUse oldSupportSolveUse = SUPPORT_SOLVE_UNSPECIFIED;
651
652 belosOp->uninitialize( &oldLP, NULL, &oldIterSolver, &oldFwdOpSrc,
653 NULL, &oldIsExternalPrec, &oldApproxFwdOpSrc, &oldSupportSolveUse );
654
655 //
656 // Create the Belos linear problem
657 // NOTE: If one exists already, reuse it.
658 //
659
661 RCP<LP_t> lp;
662 if (oldLP != Teuchos::null) {
663 lp = oldLP;
664 }
665 else {
666 lp = rcp(new LP_t());
667 }
668
669 //
670 // Set the operator
671 //
672
673 lp->setOperator(fwdOp);
674
675 //
676 // Set the preconditioner
677 //
678
679 if(prec.get()) {
680 RCP<const LinearOpBase<Scalar> > unspecified = prec->getUnspecifiedPrecOp();
681 RCP<const LinearOpBase<Scalar> > left = prec->getLeftPrecOp();
682 RCP<const LinearOpBase<Scalar> > right = prec->getRightPrecOp();
683 TEUCHOS_TEST_FOR_EXCEPTION(
684 !( left.get() || right.get() || unspecified.get() ), std::logic_error
685 ,"Error, at least one preconditoner linear operator objects must be set!"
686 );
687 if (nonnull(unspecified)) {
688 if (paramList_->get<bool>(LeftPreconditionerIfUnspecified_name, false))
689 lp->setLeftPrec(unspecified);
690 else
691 lp->setRightPrec(unspecified);
692 }
693 else if (nonnull(left)) {
694 lp->setLeftPrec(left);
695 }
696 else if (nonnull(right)) {
697 lp->setRightPrec(right);
698 }
699 else {
700 // Set a left, right or split preconditioner
701 TEUCHOS_TEST_FOR_EXCEPTION(
702 nonnull(left) && nonnull(right),std::logic_error
703 ,"Error, we can not currently handle split preconditioners!"
704 );
705 }
706 }
707 if(myPrec.get()) {
708 set_extra_data<RCP<PreconditionerBase<Scalar> > >(myPrec,"Belos::InternalPrec",
709 Teuchos::inOutArg(lp), Teuchos::POST_DESTROY, false);
710 }
711 else if(prec.get()) {
712 set_extra_data<RCP<const PreconditionerBase<Scalar> > >(prec,"Belos::ExternalPrec",
713 Teuchos::inOutArg(lp), Teuchos::POST_DESTROY, false);
714 }
715
716 //
717 // Generate the parameter list.
718 //
719
720 typedef Belos::SolverManager<Scalar,MV_t,LO_t> IterativeSolver_t;
721 RCP<IterativeSolver_t> iterativeSolver = Teuchos::null;
722 RCP<Teuchos::ParameterList> solverPL = Teuchos::rcp( new Teuchos::ParameterList() );
723
724 switch(solverType_) {
725 case SOLVER_TYPE_BLOCK_GMRES:
726 {
727 // Set the PL
728 if(paramList_.get()) {
729 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
730 Teuchos::ParameterList &gmresPL = solverTypesPL.sublist(BlockGMRES_name);
731 solverPL = Teuchos::rcp( &gmresPL, false );
732 }
733 // Create the solver
734 if (oldIterSolver != Teuchos::null) {
735 iterativeSolver = oldIterSolver;
736 iterativeSolver->setProblem( lp );
737 iterativeSolver->setParameters( solverPL );
738 }
739 else {
740 iterativeSolver = rcp(new Belos::BlockGmresSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
741 }
742 break;
743 }
744 case SOLVER_TYPE_PSEUDO_BLOCK_GMRES:
745 {
746 // Set the PL
747 if(paramList_.get()) {
748 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
749 Teuchos::ParameterList &pbgmresPL = solverTypesPL.sublist(PseudoBlockGMRES_name);
750 solverPL = Teuchos::rcp( &pbgmresPL, false );
751 }
752 //
753 // Create the solver
754 //
755 if (oldIterSolver != Teuchos::null) {
756 iterativeSolver = oldIterSolver;
757 iterativeSolver->setProblem( lp );
758 iterativeSolver->setParameters( solverPL );
759 }
760 else {
761 iterativeSolver = rcp(new Belos::PseudoBlockGmresSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
762 }
763 break;
764 }
765 case SOLVER_TYPE_BLOCK_CG:
766 {
767 // Set the PL
768 if(paramList_.get()) {
769 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
770 Teuchos::ParameterList &cgPL = solverTypesPL.sublist(BlockCG_name);
771 solverPL = Teuchos::rcp( &cgPL, false );
772 }
773 // Create the solver
774 if (oldIterSolver != Teuchos::null) {
775 iterativeSolver = oldIterSolver;
776 iterativeSolver->setProblem( lp );
777 iterativeSolver->setParameters( solverPL );
778 }
779 else {
780 iterativeSolver = rcp(new Belos::BlockCGSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
781 }
782 break;
783 }
784 case SOLVER_TYPE_PSEUDO_BLOCK_CG:
785 {
786 // Set the PL
787 if(paramList_.get()) {
788 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
789 Teuchos::ParameterList &pbcgPL = solverTypesPL.sublist(PseudoBlockCG_name);
790 solverPL = Teuchos::rcp( &pbcgPL, false );
791 }
792 //
793 // Create the solver
794 //
795 if (oldIterSolver != Teuchos::null) {
796 iterativeSolver = oldIterSolver;
797 iterativeSolver->setProblem( lp );
798 iterativeSolver->setParameters( solverPL );
799 }
800 else {
801 iterativeSolver = rcp(new Belos::PseudoBlockCGSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
802 }
803 break;
804 }
805 case SOLVER_TYPE_PSEUDO_BLOCK_STOCHASTIC_CG:
806 {
807 // Set the PL
808 if(paramList_.get()) {
809 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
810 Teuchos::ParameterList &pbcgPL = solverTypesPL.sublist(PseudoBlockStochasticCG_name);
811 solverPL = Teuchos::rcp( &pbcgPL, false );
812 }
813 //
814 // Create the solver
815 //
816 if (oldIterSolver != Teuchos::null) {
817 iterativeSolver = oldIterSolver;
818 iterativeSolver->setProblem( lp );
819 iterativeSolver->setParameters( solverPL );
820 }
821 else {
822 iterativeSolver = rcp(new Belos::PseudoBlockStochasticCGSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
823 }
824 break;
825 }
826 case SOLVER_TYPE_GCRODR:
827 {
828 // Set the PL
829 if(paramList_.get()) {
830 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
831 Teuchos::ParameterList &gcrodrPL = solverTypesPL.sublist(GCRODR_name);
832 solverPL = Teuchos::rcp( &gcrodrPL, false );
833 }
834 // Create the solver
835 if (oldIterSolver != Teuchos::null) {
836 iterativeSolver = oldIterSolver;
837 iterativeSolver->setProblem( lp );
838 iterativeSolver->setParameters( solverPL );
839 }
840 else {
841 iterativeSolver = rcp(new Belos::GCRODRSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
842 }
843 break;
844 }
845 case SOLVER_TYPE_RCG:
846 {
847 // Set the PL
848 if(paramList_.get()) {
849 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
850 Teuchos::ParameterList &rcgPL = solverTypesPL.sublist(RCG_name);
851 solverPL = Teuchos::rcp( &rcgPL, false );
852 }
853 // Create the solver
854 if (oldIterSolver != Teuchos::null) {
855 iterativeSolver = oldIterSolver;
856 iterativeSolver->setProblem( lp );
857 iterativeSolver->setParameters( solverPL );
858 }
859 else {
860 iterativeSolver = rcp(new Belos::RCGSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
861 }
862 break;
863 }
864 case SOLVER_TYPE_MINRES:
865 {
866 // Set the PL
867 if(paramList_.get()) {
868 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
869 Teuchos::ParameterList &minresPL = solverTypesPL.sublist(MINRES_name);
870 solverPL = Teuchos::rcp( &minresPL, false );
871 }
872 // Create the solver
873 if (oldIterSolver != Teuchos::null) {
874 iterativeSolver = oldIterSolver;
875 iterativeSolver->setProblem( lp );
876 iterativeSolver->setParameters( solverPL );
877 }
878 else {
879 iterativeSolver = rcp(new Belos::MinresSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
880 }
881 break;
882 }
883 case SOLVER_TYPE_TFQMR:
884 {
885 // Set the PL
886 if(paramList_.get()) {
887 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
888 Teuchos::ParameterList &minresPL = solverTypesPL.sublist(TFQMR_name);
889 solverPL = Teuchos::rcp( &minresPL, 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::TFQMRSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
899 }
900 break;
901 }
902 case SOLVER_TYPE_BICGSTAB:
903 {
904 // Set the PL
905 if(paramList_.get()) {
906 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
907 Teuchos::ParameterList &bicgstabPL = solverTypesPL.sublist(BiCGStab_name);
908 solverPL = Teuchos::rcp( &bicgstabPL, 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::BiCGStabSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
918 }
919 break;
920 }
921 case SOLVER_TYPE_FIXEDPOINT:
922 {
923 // Set the PL
924 if(paramList_.get()) {
925 Teuchos::ParameterList &solverTypesPL = paramList_->sublist(SolverTypes_name);
926 Teuchos::ParameterList &fixedPointPL = solverTypesPL.sublist(FixedPoint_name);
927 solverPL = Teuchos::rcp( &fixedPointPL, 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::FixedPointSolMgr<Scalar,MV_t,LO_t>(lp,solverPL));
937 }
938 break;
939 }
940
941 default:
942 {
943 TEUCHOS_TEST_FOR_EXCEPT(true);
944 }
945 }
946
947 //
948 // Initialize the LOWS object
949 //
950
951 belosOp->initialize(
952 lp, solverPL, iterativeSolver,
953 fwdOpSrc, prec, myPrec.get()==NULL, approxFwdOpSrc,
954 supportSolveUse, convergenceTestFrequency_
955 );
956 belosOp->setOStream(out);
957 belosOp->setVerbLevel(verbLevel);
958#ifdef TEUCHOS_DEBUG
959 if(paramList_.get()) {
960 // Make sure we read the list correctly
961 paramList_->validateParameters(*this->getValidParameters(),1); // Validate 0th and 1st level deep
962 }
963#endif
964 if(out.get() && static_cast<int>(verbLevel) > static_cast<int>(Teuchos::VERB_LOW))
965 *out << "\nLeaving Thyra::BelosLinearOpWithSolveFactory<"<<ST::name()<<">::initializeOpImpl(...) ...\n";
966
967}
968
969
970} // namespace Thyra
971
972
973#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
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
LinearOpWithSolveFactoryBase subclass implemented in terms of Belos.
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.
RCP< const LinearOpSourceBase< Scalar > > extract_approxFwdOpSrc()
RCP< const PreconditionerBase< Scalar > > extract_prec()
RCP< const LinearOpSourceBase< Scalar > > extract_fwdOpSrc()

Generated on Fri Mar 10 2023 07:13:22 for Stratimikos by doxygen 1.9.4