45 #ifndef ROL_FLETCHER_H
46 #define ROL_FLETCHER_H
53 #include "ROL_Ptr.hpp"
101 Ptr<Vector<Real> >
w_;
102 Ptr<Vector<Real> >
v_;
130 const Ptr<Constraint<Real> >
con_;
131 const Ptr<const Vector<Real> >
x_;
142 con_->applyAdjointJacobian(*(Hvp.
get(0)), *(vp.
get(1)), *
x_, tol);
143 Hvp.
get(0)->plus(*(vp.
get(0)));
145 con_->applyJacobian(*(Hvp.
get(1)), *(vp.
get(0)), *
x_, tol);
152 const Ptr<Constraint<Real> >
con_;
153 const Ptr<const Vector<Real> >
x_;
168 con_->applyPreconditioner(*(Hvp.
get(1)),*(vp.
get(1)),*
x_,*
x_, zero);
177 ROL::ParameterList &parlist)
181 y_ = conVec.
dual().clone();
182 g_ = optVec.
dual().clone();
188 w_ = optVec.
dual().clone();
189 v_ = conVec.
dual().clone();
200 vv_ = makePtr<PartitionedVector<Real>>(std::vector<Ptr<Vector<Real>> >({
v1_,
v2_}));
204 ww_ = makePtr<PartitionedVector<Real>>(std::vector<Ptr<Vector<Real>> >({
w1_,
w2_}));
208 bb_ = makePtr<PartitionedVector<Real>>(std::vector<Ptr<Vector<Real>> >({
b1_,
b2_}));
210 ROL::ParameterList& sublist = parlist.sublist(
"Step").sublist(
"Fletcher");
211 HessianApprox_ = sublist.get(
"Level of Hessian Approximation", 0);
215 delta_ = sublist.get(
"Regularization Parameter", 0.0);
217 useInexact_ = sublist.get(
"Inexact Solves",
false);
219 ROL::ParameterList krylovList;
220 Real atol =
static_cast<Real
>(1e-12);
221 Real rtol =
static_cast<Real
>(1e-2);
222 krylovList.sublist(
"General").sublist(
"Krylov").set(
"Type",
"GMRES");
223 krylovList.sublist(
"General").sublist(
"Krylov").set(
"Absolute Tolerance", atol);
224 krylovList.sublist(
"General").sublist(
"Krylov").set(
"Relative Tolerance", rtol);
225 krylovList.sublist(
"General").sublist(
"Krylov").set(
"Iteration Limit", 200);
226 krylov_ = KrylovFactory<Real>(krylovList);
230 obj_->update(x,flag,iter);
231 con_->update(x,flag,iter);
291 con_->applyAdjointHessian( *
gPhi_, *
y_, *
wg_, x, tol2 ); tol2 = origTol;
294 obj_->hessVec( *
Tv_, *
wg_, x, tol2 ); tol2 = origTol;
295 gPhi_->axpy( static_cast<Real>(-1), *
Tv_ );
297 con_->applyAdjointHessian( *
Tv_, *
vg_, *
gL_, x, tol2 ); tol2 = origTol;
303 con_->applyAdjointJacobian( *
Tv_, *
c_, x, tol2 ); tol2 = origTol;
327 obj_->hessVec( hv, v, x, tol2 ); tol2 = origTol;
328 con_->applyAdjointHessian( *
Tv_, *
y_, v, x, tol2 ); tol2 = origTol;
329 hv.
axpy(static_cast<Real>(-1), *
Tv_ );
333 hv.
scale( static_cast<Real>(-1) );
341 obj_->hessVec( *
Tv_, *
w_, x, tol2 ); tol2 = origTol;
343 con_->applyAdjointHessian( *
Tv_, *
y_, *
w_, x, tol2 ); tol2 = origTol;
344 hv.
axpy( static_cast<Real>(-1), *
Tv_ );
349 con_->applyJacobian( *
b2_, v, x, tol2 ); tol2 = origTol;
350 con_->applyAdjointJacobian( *
Tv_, *
b2_, x, tol2 ); tol2 = origTol;
352 con_->applyAdjointHessian( *
Tv_, *
c_, v, x, tol2); tol2 = origTol;
364 bool refine =
false) {
366 ROL::Ptr<LinearOperator<Real> > K
367 = ROL::makePtr<AugSystem>(
con_, makePtrFromRef(x),
delta_);
368 ROL::Ptr<LinearOperator<Real> > P
369 = ROL::makePtr<AugSystemPrecond>(
con_, makePtrFromRef(x));
379 K->apply(*
vv_, *
ww_, tol); tol = origTol;
381 b1_->axpy( static_cast<Real>(-1), *
v1_ );
382 b2_->axpy( static_cast<Real>(-1), *
v2_ );
390 krylov_->resetAbsoluteTolerance(tol);
Provides the interface to evaluate objective functions.
Real value(const Vector< Real > &x, Real &tol)
Compute value.
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
Ptr< Vector< Real > > w2_
Ptr< Vector< Real > > gPhi_
virtual void scale(const Real alpha)=0
Compute where .
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
const Ptr< Constraint< Real > > con_
Real quadPenaltyParameter_
ROL::Ptr< const Vector< Real > > get(size_type i) const
virtual void plus(const Vector &x)=0
Compute , where .
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
Defines the linear algebra of vector space on a generic partitioned vector.
void objValue(const Vector< Real > &x, Real &tol)
const Ptr< Objective< Real > > obj_
void applyInverse(Vector< Real > &Hv, const Vector< Real > &v, Real &tol) const
Apply inverse of linear operator.
Ptr< Vector< Real > > vg_
Contains definitions of custom data types in ROL.
void solveAugmentedSystem(Vector< Real > &v1, Vector< Real > &v2, const Vector< Real > &b1, const Vector< Real > &b2, const Vector< Real > &x, Real &tol, bool refine=false)
Fletcher(const ROL::Ptr< Objective< Real > > &obj, const ROL::Ptr< Constraint< Real > > &con, const Vector< Real > &optVec, const Vector< Real > &conVec, ROL::ParameterList &parlist)
Ptr< Vector< Real > > wg_
const Ptr< const Vector< Real > > x_
const Ptr< const Vector< Real > > x_
AugSystemPrecond(const Ptr< Constraint< Real > > con, const Ptr< const Vector< Real > > x)
Ptr< Vector< Real > > czeros_
Ptr< Vector< Real > > w1_
Ptr< Vector< Real > > Tv_
Defines the linear algebra or vector space interface.
Ptr< Krylov< Real > > krylov_
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
Ptr< Vector< Real > > gL_
Ptr< Vector< Real > > v1_
Ptr< Vector< Real > > b1_
AugSystem(const Ptr< Constraint< Real > > &con, const Ptr< const Vector< Real > > &x, const Real delta)
bool isMultiplierComputed_
Ptr< Vector< Real > > scaledc_
void apply(Vector< Real > &Hv, const Vector< Real > &v, Real &tol) const
Apply linear operator.
void conValue(const Vector< Real > &x, Real &tol)
void set(const V &x)
Set where .
void computeMultipliers(const Vector< Real > &x, const Real tol)
Ptr< Vector< Real > > xzeros_
Ptr< PartitionedVector< Real > > vv_
Provides the interface to apply a linear operator.
const Ptr< Constraint< Real > > con_
void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
virtual void set(const Vector &x)
Set where .
void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update objective function.
Ptr< Vector< Real > > v2_
Ptr< Vector< Real > > b2_
Ptr< PartitionedVector< Real > > bb_
const Ptr< Constraint< Real > > con_
Ptr< PartitionedVector< Real > > ww_
void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
void objGrad(const Vector< Real > &x, Real &tol)
Defines the general constraint operator interface.
void apply(Vector< Real > &Hv, const Vector< Real > &v, Real &tol) const
Apply linear operator.