Solves a system of linear equations with general minimal residuals (GMRes).
More...
#include <gmres.hh>
|
virtual const uint | dimX () const |
| Returns the size of the image space of the operator (number of rows of the corresponding matrix) More...
|
|
virtual const uint | dimY () const |
| Returns the size of the source space of the operator (number of columns of the corresponding matrix) More...
|
|
Real | epsilon () const |
| Returns the residual. More...
|
|
| GMRes (Operator< F > &A, Real maxeps, int maxit=0, uint rs=0, uint relres=0, Operator< F > *W=0) |
| Constructor. More...
|
|
uint | iterations () const |
| Returns the number of iterations. More...
|
|
void | operator() () |
| Application method without second argument. Used for parallel solvers. More...
|
|
virtual void | operator() (const Function< c_type > &fncY, Function< c_type > &fncX) |
| Application operator for complex function fncY . More...
|
|
virtual void | operator() (const Function< r_type > &fncY, Function< F > &fncX) |
| Application operator for real function fncY . More...
|
|
void | operator() (const Matrix< c_type > &mX, Matrix< c_type > &mY) |
| Application method to complex matrices. Calls apply_() More...
|
|
void | operator() (const Matrix< r_type > &mX, Matrix< F > &mY) |
| Application method to real matrices. Calls function apply() More...
|
|
virtual void | operator() (const Vector< c_type > &fncY, Vector< c_type > &fncX) |
| Application operator for complex function fncY . More...
|
|
virtual void | operator() (const Vector< r_type > &fncY, Vector< F > &fncX) |
| Application operator for real vector fncY . More...
|
|
virtual void | show_messages () |
|
virtual | ~GMRes () |
|
|
virtual std::ostream & | info (std::ostream &os) const |
|
|
virtual void | apply_ () throw (NoConvergence, MissingFeature) |
| Intrinsic application method without argument. More...
|
|
virtual void | apply_ (const Vector< F > &fncY, Vector< F > &fncX) throw (NoConvergence, MissingFeature) |
| Intrinsic application method, i.e. More...
|
|
template<class F>
class concepts::GMRes< F >
Solves a system of linear equations with general minimal residuals (GMRes).
- Examples
- linearDG1d.cc.
Definition at line 24 of file gmres.hh.
◆ c_type
◆ r_type
◆ type
◆ GMRes()
Constructor.
- Parameters
-
A | A matrix |
maxeps | Maximal residual |
maxit | Maximal number of iterations |
rs | Restart |
relres | 0: absolute residual, 1: relative residual |
W | Optional preconditioner |
Definition at line 34 of file gmres.hh.
◆ ~GMRes()
◆ apply_() [1/2]
◆ apply_() [2/2]
◆ dimX()
Returns the size of the image space of the operator (number of rows of the corresponding matrix)
Definition at line 93 of file compositions.hh.
◆ dimY()
Returns the size of the source space of the operator (number of columns of the corresponding matrix)
Definition at line 98 of file compositions.hh.
◆ epsilon()
Returns the residual.
Calling this method makes only sence after a linear system has been solved.
Definition at line 48 of file gmres.hh.
◆ info()
template<class F >
virtual std::ostream& concepts::GMRes< F >::info |
( |
std::ostream & |
os | ) |
const |
|
protectedvirtual |
◆ iterations()
Returns the number of iterations.
Calling this method makes only sence after a linear system has been solved.
Definition at line 43 of file gmres.hh.
◆ operator()() [1/7]
◆ operator()() [2/7]
Application operator for complex function fncY
.
Computes fncX
= A(fncY
) where A is this operator. fncX
becomes complex.
In derived classes its enough to implement the operator() for complex Operator's. If a real counterpart is not implemented, the function fncY
is splitted into real and imaginary part and the application operator for real functions is called for each. Then the result is combined.
If in a derived class the operator() for complex Operator's is not implemented, a exception is thrown from here.
Reimplemented from concepts::Operator< F >.
◆ operator()() [3/7]
Application operator for real function fncY
.
Computes fncX
= A(fncY
) where A is this operator.
fncX
becomes the type of the operator, for real data it becomes real, for complex data it becomes complex.
In derived classes its enough to implement the operator() for real Operator's. If a complex counterpart is not implemented, the function fncY
is transformed to a complex function and then the application operator for complex functions is called.
If in a derived class the operator() for real Operator's is not implemented, a exception is thrown from here.
Reimplemented from concepts::Operator< F >.
◆ operator()() [4/7]
Application method to complex matrices. Calls apply_()
◆ operator()() [5/7]
Application method to real matrices. Calls function apply()
◆ operator()() [6/7]
Application operator for complex function fncY
.
Computes fncX
= A(fncY
) where A is this operator. fncX
becomes complex.
In derived classes its enough to implement the operator() for complex Operator's. If a real counterpart is not implemented, the vector fncY
is splitted into real and imaginary part and the application operator for real vectors is called for each. Then the result is combined
If in a derived class the operator() for complex Operator's i not implemented, a exception is thrown from here.
◆ operator()() [7/7]
Application operator for real vector fncY
.
Computes fncX
= A(fncY
) where A is this operator.
Type of fncX
becomes that of the operator, for real data it becomes real, for complex data it becomes complex.
In derived classes its enough to implement the operator() for real Operator's. If a complex counterpart is not implemented, the vector fncY
is transformed to a complex vector and then the application for complex vectors is called.
If in a derived class the operator() for real Operator's is not implemented, a exception is thrown from here.
◆ show_messages()
◆ A_
◆ dimX_
Dimension of image space and the source space.
Definition at line 104 of file compositions.hh.
◆ dimY_
◆ eps_
Current residual.
Definition at line 78 of file gmres.hh.
◆ it_
Number of iterations.
Definition at line 80 of file gmres.hh.
◆ maxeps_
Convergence criterion.
Definition at line 74 of file gmres.hh.
◆ maxit_
Maximal number of iterations until abortion.
Definition at line 76 of file gmres.hh.
◆ relres_
relres = 0: absolute residual, relres = 1: relative residual
Definition at line 84 of file gmres.hh.
◆ rs_
◆ W_
Optional preconditioner.
Definition at line 69 of file gmres.hh.
The documentation for this class was generated from the following file: