Belos Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
BelosGmresIteration.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_GMRES_ITERATION_HPP
11 #define BELOS_GMRES_ITERATION_HPP
12 
17 #include "BelosConfigDefs.hpp"
18 #include "BelosTypes.hpp"
19 #include "BelosIteration.hpp"
20 
21 namespace Belos {
22 
24 
25 
30  template <class ScalarType, class MV>
36  int curDim;
37 
40 
43 
57 
58 
59  GmresIterationState() : curDim(0), V(Teuchos::null), Z(Teuchos::null),
60  H(Teuchos::null), R(Teuchos::null),
61  z(Teuchos::null)
62  {}
63  };
64 
66 
68 
69 
81  class GmresIterationInitFailure : public BelosError {public:
82  GmresIterationInitFailure(const std::string& what_arg) : BelosError(what_arg)
83  {}};
84 
91  class GmresIterationOrthoFailure : public BelosError {public:
92  GmresIterationOrthoFailure(const std::string& what_arg) : BelosError(what_arg)
93  {}};
94 
101  class GmresIterationLAPACKFailure : public BelosError {public:
102  GmresIterationLAPACKFailure(const std::string& what_arg) : BelosError(what_arg)
103  {}};
104 
106 
107 
108 template<class ScalarType, class MV, class OP>
109 class GmresIteration : virtual public Iteration<ScalarType,MV,OP> {
110 
111  public:
112 
114 
115 
129  virtual void initializeGmres(GmresIterationState<ScalarType,MV>& newstate) = 0;
130 
137  virtual GmresIterationState<ScalarType,MV> getState() const = 0;
139 
141 
142 
144 
147  virtual void updateLSQR( int dim = -1 ) = 0;
148 
150  virtual int getCurSubspaceDim() const = 0;
151 
153  virtual int getMaxSubspaceDim() const = 0;
154 
156 
158 
159 
166  virtual void setSize(int blockSize, int numBlocks) = 0;
167 
169 
170 };
171 
172 } // end Belos namespace
173 
174 #endif /* BELOS_GMRES_ITERATION_HPP */
Collection of types and exceptions used within the Belos solvers.
virtual void updateLSQR(int dim=-1)=0
Method for updating QR factorization of upper Hessenberg matrix.
Teuchos::RCP< const MV > V
The current Krylov basis.
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > z
The current right-hand side of the least squares system RY = Z.
GmresIterationLAPACKFailure is thrown when a nonzero return value is passed back from an LAPACK routi...
virtual int getCurSubspaceDim() const =0
Get the dimension of the search subspace used to generate the current solution to the linear problem...
Structure to contain pointers to GmresIteration state variables.
Pure virtual base class which describes the basic interface to the linear solver iteration.
GmresIterationInitFailure is thrown when the GmresIteration object is unable to generate an initial i...
Teuchos::RCP< const MV > Z
The current preconditioned Krylov basis (only used in flexible GMRES).
virtual int getMaxSubspaceDim() const =0
Get the maximum dimension allocated for the search subspace.
int curDim
The current dimension of the reduction.
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > R
The current upper-triangular matrix from the QR reduction of H.
virtual GmresIterationState< ScalarType, MV > getState() const =0
Get the current state of the linear solver.
virtual void initializeGmres(GmresIterationState< ScalarType, MV > &newstate)=0
Initialize the solver to an iterate, providing a complete state.
GmresIterationLAPACKFailure(const std::string &what_arg)
GmresIterationInitFailure(const std::string &what_arg)
Parent class to all Belos exceptions.
Definition: BelosTypes.hpp:28
GmresIterationOrthoFailure(const std::string &what_arg)
GmresIterationOrthoFailure is thrown when the GmresIteration object is unable to compute independent ...
Belos header file which uses auto-configuration information to include necessary C++ headers...
virtual void setSize(int blockSize, int numBlocks)=0
Set the blocksize and number of blocks to be used by the iterative solver in solving this linear prob...
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > H
The current Hessenberg matrix.