ROL
ROL_ConjugateResiduals.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ************************************************************************
3 //
4 // Rapid Optimization Library (ROL) Package
5 // Copyright (2014) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact lead developers:
38 // Drew Kouri (dpkouri@sandia.gov) and
39 // Denis Ridzal (dridzal@sandia.gov)
40 //
41 // ************************************************************************
42 // @HEADER
43 
44 #ifndef ROL_CONJUGATERESIDUALS_H
45 #define ROL_CONJUGATERESIDUALS_H
46 
51 #include "ROL_Krylov.hpp"
52 #include "ROL_LinearOperator.hpp"
53 #include "ROL_Vector.hpp"
54 #include "ROL_Types.hpp"
55 
56 namespace ROL {
57 
58 template<class Real>
59 class ConjugateResiduals : public Krylov<Real> {
60 
63  ROL::Ptr<Vector<Real> > r_;
64  ROL::Ptr<Vector<Real> > v_;
65  ROL::Ptr<Vector<Real> > p_;
66  ROL::Ptr<Vector<Real> > Ap_;
67  ROL::Ptr<Vector<Real> > MAp_;
68 
69 public:
70  ConjugateResiduals( Real absTol = 1.e-4, Real relTol = 1.e-2, int maxit = 100, bool useInexact = false )
71  : Krylov<Real>(absTol,relTol,maxit), isInitialized_(false), useInexact_(useInexact) {}
72 
73  // Run Krylov Method
75  int &iter, int &flag ) {
76  if ( !isInitialized_ ) {
77  r_ = x.clone();
78  v_ = b.clone();
79  p_ = x.clone();
80  Ap_ = b.clone();
81  MAp_ = x.clone();
82  isInitialized_ = true;
83  }
84 
85  // Initialize
86  Real rnorm = b.norm();
88  Real itol = std::sqrt(ROL_EPSILON<Real>());
89  x.zero();
90 
91  // Apply preconditioner to residual
92  M.applyInverse(*r_,b,itol);
93 
94  // Initialize direction p
95  p_->set(*r_);
96 
97  // Get Hessian tolerance
98  if ( useInexact_ ) {
99  itol = rtol/((Real)Krylov<Real>::getMaximumIteration() * rnorm);
100  }
101 
102  // Apply Hessian to residual
103  A.apply(*v_, *r_, itol);
104 
105  // Apply Hessian to direction p
106  //A.apply(*Ap_, *p_, itol);
107  Ap_->set(*v_);
108 
109  // Initialize scalar quantities
110  iter = 0;
111  flag = 0;
112  Real kappa(0), beta(0), alpha(0), tmp(0);
113  Real gHg = r_->dot(v_->dual());
114 
115  for (iter = 0; iter < (int)Krylov<Real>::getMaximumIteration(); iter++) {
116  itol = std::sqrt(ROL_EPSILON<Real>());
117  M.applyInverse(*MAp_, *Ap_, itol);
118  kappa = MAp_->dot(Ap_->dual());
119  //if ( gHg <= 0.0 || kappa <= 0.0 ) {
120  //flag = 2;
121  //break;
122  //}
123  alpha = gHg/kappa;
124 
125  x.axpy(alpha,*p_);
126 
127  r_->axpy(-alpha,*MAp_);
128  rnorm = r_->norm();
129  if ( rnorm < rtol ) {
130  break;
131  }
132 
133  if ( useInexact_ ) {
134  itol = rtol/((Real)Krylov<Real>::getMaximumIteration() * rnorm);
135  }
136  A.apply(*v_, *r_, itol);
137  tmp = gHg;
138  gHg = r_->dot(v_->dual());
139  beta = gHg/tmp;
140 
141  p_->scale(beta);
142  p_->plus(*r_);
143 
144  Ap_->scale(beta);
145  Ap_->plus(*v_);
146  }
147  if ( iter == (int)Krylov<Real>::getMaximumIteration() ) {
148  flag = 1;
149  }
150  else {
151  iter++;
152  }
153  return rnorm;
154  }
155 };
156 
157 
158 }
159 
160 #endif
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
Provides definition of the Conjugate Residual solver.
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
Definition: ROL_Vector.hpp:153
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.
virtual void zero()
Set to zero vector.
Definition: ROL_Vector.hpp:167
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:80
ROL::Ptr< Vector< Real > > r_
ROL::Ptr< Vector< Real > > Ap_
virtual void applyInverse(Vector< Real > &Hv, const Vector< Real > &v, Real &tol) const
Apply inverse of linear operator.
Real run(Vector< Real > &x, LinearOperator< Real > &A, const Vector< Real > &b, LinearOperator< Real > &M, int &iter, int &flag)
Provides definitions for Krylov solvers.
Definition: ROL_Krylov.hpp:58
Provides the interface to apply a linear operator.
ROL::Ptr< Vector< Real > > v_
ConjugateResiduals(Real absTol=1.e-4, Real relTol=1.e-2, int maxit=100, bool useInexact=false)
ROL::Ptr< Vector< Real > > p_
ROL::Ptr< Vector< Real > > MAp_
virtual Real norm() const =0
Returns where .