Belos Version of the Day
Loading...
Searching...
No Matches
BelosInnerSolver.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Belos: Block Linear Solvers Package
4//
5// Copyright 2004-2016 NTESS and the Belos contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef __Belos_InnerSolver_hpp
11#define __Belos_InnerSolver_hpp
12
17#include <Teuchos_ScalarTraits.hpp>
18
19namespace Belos {
20
32 template<class Scalar, class MV, class OP>
33 Teuchos::RCP<LinearProblem<Scalar, MV, OP> >
35 const Teuchos::RCP<const MV>& B)
36 {
37 using Teuchos::is_null;
38 using Teuchos::nonnull;
39 using Teuchos::null;
40 using Teuchos::RCP;
41 using Teuchos::rcp;
44
45 RCP<const OP> A = problem->getOperator ();
46 RCP<MV> X_orig = problem->getLHS ();
47 TEUCHOS_TEST_FOR_EXCEPTION(is_null (X_orig), std::invalid_argument,
48 "problemWithNewRHS(): The original LinearProblem's "
49 "initial guess / current approximate solution (getLHS())"
50 " is null. We need an initial guess or current approxim"
51 "ate solution in order to know the domain of the (right-"
52 "preconditioned, if applicable) operator. This is "
53 "because Belos::MultiVecTraits does not include the idea"
54 " of the domain and range of an operator, or the space "
55 "to which a vector belongs.");
56 TEUCHOS_TEST_FOR_EXCEPTION(is_null (B), std::invalid_argument,
57 "problemWithNewRHS(): the given new right-hand side B "
58 "is null.");
59 RCP<MV> X = MVT::CloneCopy (problem->getLHS ());
60
61 RCP<lp_type> lp (new lp_type (A, X, B));
62 lp->setLeftPrec (problem->getLeftPrec ());
63 lp->setRightPrec (problem->getRightPrec ());
64 // Compute initial residual(s) and prepare the problem for solution.
65 lp->setProblem ();
66 return lp;
67 }
68
69
105 template<class Scalar, class MV, class OP>
107 public:
109 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType magnitude_type;
112
114 virtual ~InnerSolver() {}
115
122 virtual Teuchos::RCP<const Teuchos::ParameterList>
124
168 virtual InnerSolveResult
169 solve (const Teuchos::RCP<MV>& X,
170 const Teuchos::RCP<const MV>& B,
171 const magnitude_type convTol,
172 const int maxItersPerRestart,
173 const int maxNumRestarts) = 0;
174
218 virtual InnerSolveResult
219 solve (const Teuchos::RCP<MV>& X,
220 const Teuchos::RCP<const MV>& B) = 0;
221 };
222
223
232 template<class Scalar, class MV, class OP>
234 public:
235 static void
237 const MV& x,
238 MV& y,
240 {
241 using Teuchos::RCP;
242 using Teuchos::rcpFromRef;
243
244 TEUCHOS_TEST_FOR_EXCEPTION(trans != NOTRANS, std::invalid_argument,
245 "Belos::InnerSolver is not able to solve the "
246 "transposed system.");
249 (void) Op.solve (y_ptr, x_ptr);
250 }
251
252 };
253
261 template<class Scalar, class MV, class OP>
263 public:
265 (void) Scalar::this_specialization_is_not_defined();
266 (void) MV::this_specialization_is_not_defined();
267 (void) OP::this_specialization_is_not_defined();
268 }
269 };
270
290 template<class Scalar, class MV, class OP>
292 public:
298
303 static Teuchos::RCP<OP>
305 {
306 using Teuchos::rcp;
307 using Teuchos::rcp_implicit_cast;
308 // If this class is not specialized for the given combination of
309 // (Scalar, MV, OP), the constructor of wrapper_type here will
310 // (deliberately) raise a compiler error.
312 }
313
324 static Teuchos::RCP<inner_solver_type>
325 getInnerSolver (const Teuchos::RCP<operator_type>& op)
326 {
327 using Teuchos::RCP;
328 using Teuchos::rcp_dynamic_cast;
329 // If this class is not specialized for the given combination of
330 // (Scalar, MV, OP), the instantiation of the wrapper_type class
331 // here will (deliberately) raise a compiler error.
333 return wrapper->getInnerSolver();
334 }
335 };
336
337} // namespace Belos
338
339#endif // __Belos_InnerSolver_hpp
Class which describes the linear problem to be solved by the iterative solver.
Declaration of basic traits for the multivector type.
Class which defines basic traits for the operator type.
Represents the result of an inner solve.
Wrap an InnerSolver in an OP (operator).
InnerSolver< scalar_type, multivector_type, operator_type > inner_solver_type
static Teuchos::RCP< inner_solver_type > getInnerSolver(const Teuchos::RCP< operator_type > &op)
Return the given wrapper's inner solver object.
UndefinedWrapperType< Scalar, MV, OP > wrapper_type
static Teuchos::RCP< OP > makeInnerSolverOperator(const Teuchos::RCP< InnerSolver< Scalar, MV, OP > > &solver)
Wrap the given inner solver in a wrapper_type.
Inner solver interface.
virtual ~InnerSolver()
Virtual destructor, for correctness.
virtual Teuchos::RCP< const Teuchos::ParameterList > getCurrentParameters() const =0
Current parameters for the inner solver implementation.
Teuchos::ScalarTraits< Scalar >::magnitudeType magnitude_type
virtual InnerSolveResult solve(const Teuchos::RCP< MV > &X, const Teuchos::RCP< const MV > &B)=0
Solve for the given right-hand side(s) B.
virtual InnerSolveResult solve(const Teuchos::RCP< MV > &X, const Teuchos::RCP< const MV > &B, const magnitude_type convTol, const int maxItersPerRestart, const int maxNumRestarts)=0
Solve for the given right-hand side(s) B.
static void Apply(const InnerSolver< Scalar, MV, OP > &Op, const MV &x, MV &y, ETrans trans=NOTRANS)
Class which defines basic traits for the operator type.
Alternative run-time polymorphic interface for operators.
Undefined wrapper type, to check at compile time whether InnerSolverTraits has been specialized.
ETrans
Whether to apply the (conjugate) transpose of an operator.
Teuchos::RCP< LinearProblem< Scalar, MV, OP > > problemWithNewRHS(const Teuchos::RCP< const LinearProblem< Scalar, MV, OP > > &problem, const Teuchos::RCP< const MV > &B)
New LinearProblem with different right-hand side.

Generated for Belos by doxygen 1.9.8