19#include "ROL_LAPACK.hpp"
20#include "ROL_LinearAlgebra.hpp"
28 std::vector<Ptr<Vector<Real>>>
data_;
36 for (
unsigned i = 0u; i < numVec; ++i)
45 else if (ind ==
size_) {
56 const Ptr<Vector<Real>>
get(
unsigned ind) {
86 Ptr<std::vector<Real>>
res_;
87 Ptr<VectorArray<Real>>
V_,
Z_;
105 ROL::ParameterList &gList = parlist.sublist(
"General");
106 ROL::ParameterList &kList = gList.sublist(
"Krylov");
108 useInexact_ = gList.get(
"Inexact Hessian-Times-A-Vector",
false);
112 H_ = makePtr<SDMatrix>( maxit+1, maxit );
113 cs_ = makePtr<SDVector>( maxit );
114 sn_ = makePtr<SDVector>( maxit );
115 s_ = makePtr<SDVector>( maxit+1 );
116 y_ = makePtr<SDVector>( maxit+1 );
117 cnorm_ = makePtr<SDVector>( maxit );
118 res_ = makePtr<std::vector<Real>>(maxit+1,
zero);
120 V_ = makePtr<VectorArray<Real>>();
121 Z_ = makePtr<VectorArray<Real>>();
133 Real
zero(0), one(1);
140 V_->allocate(b,std::min(20u,
static_cast<unsigned>(maxit))+1u);
141 Z_->allocate(x,std::min(20u,
static_cast<unsigned>(maxit)));
146 Real itol = std::sqrt(ROL_EPSILON<Real>());
163 (*res_)[0] =
r_->norm();
166 *
outStream_ <<
"GMRES Iteration " << 0 <<
", Residual = " << (*res_)[0] <<
"\n";
169 Real rtol = std::min(absTol,relTol*(*
res_)[0]);
170 if ((*
res_)[0] <= rtol) {
176 V_->get(0)->set(*
r_);
177 V_->get(0)->scale(one/(*
res_)[0]);
179 (*s_)(0) = (*
res_)[0];
181 for( iter=0; iter<maxit; ++iter ) {
183 if(
useInexact_ ) itol = rtol/(maxit*(*res_)[iter]);
192 for(
int k=0; k<=iter; ++k ) {
193 (*H_)(k,iter) =
w_->dot(*(
V_->get(k)));
194 w_->axpy( -(*
H_)(k,iter), *(
V_->get(k)) );
197 (*H_)(iter+1,iter) =
w_->norm();
199 (
V_->get(iter+1))->set(*
w_);
200 (
V_->get(iter+1))->scale(one/((*
H_)(iter+1,iter)));
203 for(
int k=0; k<=iter-1; ++k ) {
204 temp = (*cs_)(k)*(*
H_)(k,iter) + (*
sn_)(k)*(*
H_)(k+1,iter);
205 (*H_)(k+1,iter) = -(*
sn_)(k)*(*
H_)(k,iter) + (*
cs_)(k)*(*
H_)(k+1,iter);
206 (*H_)(k,iter) = temp;
210 if( (*
H_)(iter+1,iter) ==
zero ) {
214 else if ( std::abs((*
H_)(iter+1,iter)) > std::abs((*
H_)(iter,iter)) ) {
215 temp = (*H_)(iter,iter) / (*
H_)(iter+1,iter);
216 (*sn_)(iter) = one / std::sqrt( one + temp*temp );
217 (*cs_)(iter) = temp*(*
sn_)(iter);
220 temp = (*H_)(iter+1,iter) / (*
H_)(iter,iter);
221 (*cs_)(iter) = one / std::sqrt( one + temp*temp );
222 (*sn_)(iter) = temp*(*
cs_)(iter);
226 temp = (*cs_)(iter)*(*
s_)(iter);
227 (*s_)(iter+1) = -(*
sn_)(iter)*(*
s_)(iter);
229 (*H_)(iter,iter) = (*
cs_)(iter)*(*
H_)(iter,iter) + (*
sn_)(iter)*(*
H_)(iter+1,iter);
230 (*H_)(iter+1,iter) =
zero;
231 (*res_)[iter+1] = std::abs((*
s_)(iter+1));
234 *
outStream_ <<
"GMRES Iteration " << iter+1 <<
", Residual = " << (*res_)[iter+1] <<
"\n";
238 const char uplo =
'U';
239 const char trans =
'N';
240 const char diag =
'N';
241 const char normin =
'N';
245 lapack_.LATRS(uplo, trans, diag, normin, iter+1,
H_->values(), maxit+1,
y_->values(), &scaling,
cnorm_->values(), &info);
249 for(
int k=0; k<=iter;++k ) {
250 z_->axpy((*
y_)(k),*(
Z_->get(k)));
253 if( (*
res_)[iter+1] <= rtol ) {
264 return (*
res_)[iter];
267 return (*
res_)[iter+1];
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
Contains definitions of custom data types in ROL.
Preconditioned GMRES solver.
GMRES(ParameterList &parlist)
LA::Vector< Real > SDVector
Ptr< VectorArray< Real > > Z_
LAPACK< int, Real > lapack_
void enableOutput(std::ostream &outStream)
Ptr< std::ostream > outStream_
LA::Matrix< Real > SDMatrix
Ptr< std::vector< Real > > res_
Ptr< VectorArray< Real > > V_
Real run(Vector< Real > &x, LinearOperator< Real > &A, const Vector< Real > &b, LinearOperator< Real > &M, int &iter, int &flag)
Provides definitions for Krylov solvers.
unsigned getMaximumIteration(void) const
Real getAbsoluteTolerance(void) const
Real getRelativeTolerance(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.
std::vector< Ptr< Vector< Real > > > data_
const Ptr< Vector< Real > > get(const Vector< Real > &x, unsigned ind)
void allocate(const Vector< Real > &x, unsigned numVec)
const Ptr< Vector< Real > > get(unsigned ind)
Defines the linear algebra or vector space interface.
virtual void plus(const Vector &x)=0
Compute , where .
virtual void zero()
Set to zero vector.
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.