53 #include "ROL_LAPACK.hpp"
54 #include "ROL_LinearAlgebra.hpp"
67 ROL::Ptr<Vector<Real> >
r_;
68 ROL::Ptr<Vector<Real> >
z_;
69 ROL::Ptr<Vector<Real> >
w_;
71 ROL::Ptr<SDMatrix>
H_;
72 ROL::Ptr<SDVector>
cs_;
73 ROL::Ptr<SDVector>
sn_;
74 ROL::Ptr<SDVector>
s_;
75 ROL::Ptr<SDVector>
y_;
78 ROL::Ptr<std::vector<Real> >
res_;
96 ROL::ParameterList &gList = parlist.sublist(
"General");
97 ROL::ParameterList &kList = gList.sublist(
"Krylov");
99 useInexact_ = gList.get(
"Inexact Hessian-Times-A-Vector",
false);
103 H_ = ROL::makePtr<SDMatrix>( maxit+1, maxit );
104 cs_ = ROL::makePtr<SDVector>( maxit );
105 sn_ = ROL::makePtr<SDVector>( maxit );
106 s_ = ROL::makePtr<SDVector>( maxit+1 );
107 y_ = ROL::makePtr<SDVector>( maxit+1 );
108 cnorm_ = ROL::makePtr<SDVector>( maxit );
109 res_ = ROL::makePtr<std::vector<Real>>(maxit+1,
zero);
122 Real
zero(0), one(1);
132 Real itol = std::sqrt(ROL_EPSILON<Real>());
149 std::vector<ROL::Ptr<Vector<Real > > >
V;
150 std::vector<ROL::Ptr<Vector<Real > > > Z;
152 (*res_)[0] =
r_->norm();
155 *
outStream_ <<
"GMRES Iteration " << 0 <<
", Residual = " << (*res_)[0] <<
"\n";
159 Real rtol = std::min(absTol,relTol*(*
res_)[0]);
160 if ((*
res_)[0] <= rtol) {
166 V.push_back(b.
clone());
168 (V[0])->scale(one/(*
res_)[0]);
170 (*s_)(0) = (*
res_)[0];
172 for( iter=0; iter<maxit; ++iter ) {
175 itol = rtol/(maxit*(*res_)[iter]);
178 Z.push_back(x.
clone());
187 for(
int k=0; k<=iter; ++k ) {
188 (*H_)(k,iter) =
w_->dot(*(V[k]));
189 w_->axpy( -(*
H_)(k,iter), *(V[k]) );
192 (*H_)(iter+1,iter) =
w_->norm();
194 V.push_back( b.
clone() );
195 (V[iter+1])->set(*
w_);
196 (V[iter+1])->scale(one/((*
H_)(iter+1,iter)));
199 for(
int k=0; k<=iter-1; ++k ) {
200 temp = (*cs_)(k)*(*
H_)(k,iter) + (*
sn_)(k)*(*
H_)(k+1,iter);
201 (*H_)(k+1,iter) = -(*
sn_)(k)*(*
H_)(k,iter) + (*
cs_)(k)*(*
H_)(k+1,iter);
202 (*H_)(k,iter) = temp;
206 if( (*
H_)(iter+1,iter) ==
zero ) {
210 else if ( std::abs((*
H_)(iter+1,iter)) > std::abs((*
H_)(iter,iter)) ) {
211 temp = (*H_)(iter,iter) / (*
H_)(iter+1,iter);
212 (*sn_)(iter) = one / std::sqrt( one + temp*temp );
213 (*cs_)(iter) = temp*(*
sn_)(iter);
216 temp = (*H_)(iter+1,iter) / (*
H_)(iter,iter);
217 (*cs_)(iter) = one / std::sqrt( one + temp*temp );
218 (*sn_)(iter) = temp*(*
cs_)(iter);
222 temp = (*cs_)(iter)*(*
s_)(iter);
223 (*s_)(iter+1) = -(*
sn_)(iter)*(*
s_)(iter);
225 (*H_)(iter,iter) = (*
cs_)(iter)*(*
H_)(iter,iter) + (*
sn_)(iter)*(*
H_)(iter+1,iter);
226 (*H_)(iter+1,iter) =
zero;
227 (*res_)[iter+1] = std::abs((*
s_)(iter+1));
230 *
outStream_ <<
"GMRES Iteration " << iter+1 <<
", Residual = " << (*res_)[iter+1] <<
"\n";
234 const char uplo =
'U';
235 const char trans =
'N';
236 const char diag =
'N';
237 const char normin =
'N';
241 lapack_.LATRS(uplo, trans, diag, normin, iter+1,
H_->values(), maxit+1,
y_->values(), &scaling,
cnorm_->values(), &info);
245 for(
int k=0; k<=iter;++k ) {
246 z_->axpy((*
y_)(k),*(Z[k]));
249 if( (*res_)[iter+1] <= rtol ) {
260 return (*
res_)[iter];
263 return (*
res_)[iter+1];
277 #endif // ROL_GMRES_H
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
Real getRelativeTolerance(void) const
Real getAbsoluteTolerance(void) const
virtual void plus(const Vector &x)=0
Compute , where .
GMRES(ROL::ParameterList &parlist)
Contains definitions of custom data types in ROL.
virtual void apply(Vector< Real > &Hv, const Vector< Real > &v, Real &tol) const =0
Apply linear operator.
unsigned getMaximumIteration(void) const
virtual void zero()
Set to zero vector.
Defines the linear algebra or vector space interface.
ROL::Ptr< Vector< Real > > r_
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
ROL::Ptr< std::ostream > outStream_
LA::Matrix< Real > SDMatrix
Real run(Vector< Real > &x, LinearOperator< Real > &A, const Vector< Real > &b, LinearOperator< Real > &M, int &iter, int &flag)
ROL::Ptr< std::vector< Real > > res_
Preconditioned GMRES solver.
virtual void applyInverse(Vector< Real > &Hv, const Vector< Real > &v, Real &tol) const
Apply inverse of linear operator.
Provides definitions for Krylov solvers.
Provides the interface to apply a linear operator.
ROL::Ptr< SDVector > cnorm_
void enableOutput(std::ostream &outStream)
ROL::Ptr< Vector< Real > > w_
LA::Vector< Real > SDVector
ROL::Ptr< Vector< Real > > z_
ROL::LAPACK< int, Real > lapack_