ROL
Loading...
Searching...
No Matches
ROL_BiCGSTAB.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Rapid Optimization Library (ROL) Package
4//
5// Copyright 2014 NTESS and the ROL contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef ROL_BICGSTAB_H
11#define ROL_BICGSTAB_H
12
17#include "ROL_Krylov.hpp"
18#include "ROL_Types.hpp"
19
20namespace ROL {
21
22template<class Real>
23class BiCGSTAB : public Krylov<Real> {
24
26 const bool useInexact_;
27 Ptr<Vector<Real>> r_, r1_, p_, v_, s_, t_, h_, y_, z_;
28
29public:
30 BiCGSTAB(Real absTol = 1.e-4, Real relTol = 1.e-2, unsigned maxit = 100, bool useInexact = false)
31 : Krylov<Real>(absTol,relTol,maxit), isInitialized_(false), useInexact_(useInexact) {}
32
33 BiCGSTAB(ParameterList &parlist, bool useInexact = false)
34 : Krylov<Real>(parlist), isInitialized_(false), useInexact_(useInexact) {}
35
37 int &iter, int &flag ) {
38 if ( !isInitialized_ ) {
39 r_ = b.clone(); r1_ = b.clone(); p_ = b.clone();
40 v_ = b.clone(); s_ = b.clone(); t_ = b.clone();
41 h_ = x.clone(); y_ = x.clone(); z_ = x.clone();
42 isInitialized_ = true;
43 }
44
45 Real rho(1), rho1(1), alpha(1), beta(0), omega(1);
46 Real rnorm = b.norm();
47 Real itol = std::sqrt(ROL_EPSILON<Real>());
49 if (rnorm <= rtol) return rnorm;
50
51 x.zero();
52 v_->zero();
53 p_->zero();
54 r_->set(b);
55 r1_->set(*r_);
56
57 iter = 0;
58 flag = 0;
59
60 for (iter = 0; iter < (int)Krylov<Real>::getMaximumIteration(); iter++) {
61 rho1 = r_->dot(*r1_);
62 beta = (rho1/rho)*(alpha/omega);
63 p_->axpy(-omega,*v_);
64 p_->scale(beta);
65 p_->plus(*r_);
66
67 if ( useInexact_ )
68 itol = rtol/((Real)Krylov<Real>::getMaximumIteration() * rnorm);
69 M.applyInverse(*y_, *p_, itol);
70 A.apply(*v_, *y_, itol);
71
72 alpha = rho1 / v_->dot(*r1_);
73 h_->set(x);
74 h_->axpy(alpha,*y_);
75 s_->set(*r_);
76 s_->axpy(-alpha,*v_);
77
78 rnorm = s_->norm();
79 if (rnorm <= rtol) {
80 x.set(*h_);
81 break;
82 }
83
84 if ( useInexact_ )
85 itol = rtol/((Real)Krylov<Real>::getMaximumIteration() * rnorm);
86 M.applyInverse(*z_, *s_, itol);
87 A.apply(*t_, *z_, itol);
88
89 omega = t_->dot(*s_) / t_->dot(*t_);
90 x.set(*h_);
91 x.axpy(omega,*z_);
92 r_->set(*s_);
93 r_->axpy(-omega,*t_);
94
95 rnorm = r_->norm();
96 if (rnorm <= rtol) break;
97
98 rho = rho1;
99 }
100 if (iter == (int)Krylov<Real>::getMaximumIteration()) {
101 flag = 1;
102 }
103 else {
104 iter++;
105 }
106 return rnorm;
107 }
108};
109
110
111}
112
113#endif
Contains definitions of custom data types in ROL.
Ptr< Vector< Real > > v_
Ptr< Vector< Real > > h_
Ptr< Vector< Real > > z_
Ptr< Vector< Real > > r_
Ptr< Vector< Real > > p_
BiCGSTAB(ParameterList &parlist, bool useInexact=false)
const bool useInexact_
BiCGSTAB(Real absTol=1.e-4, Real relTol=1.e-2, unsigned maxit=100, bool useInexact=false)
Ptr< Vector< Real > > t_
Real run(Vector< Real > &x, LinearOperator< Real > &A, const Vector< Real > &b, LinearOperator< Real > &M, int &iter, int &flag)
Ptr< Vector< Real > > s_
Ptr< Vector< Real > > y_
Ptr< Vector< Real > > r1_
Provides definitions for Krylov solvers.
unsigned getMaximumIteration(void) const
Provides the interface to apply a linear operator.
virtual void apply(Vector< Real > &Hv, const Vector< Real > &v, Real &tol) const =0
Apply linear operator.
virtual void applyInverse(Vector< Real > &Hv, const Vector< Real > &v, Real &tol) const
Apply inverse of linear operator.
Defines the linear algebra or vector space interface.
virtual Real norm() const =0
Returns where .
virtual void set(const Vector &x)
Set where .
virtual void zero()
Set to zero vector.
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
virtual void axpy(const Real alpha, const Vector &x)
Compute where .