ROL
ROL_PrimalDualSystemStep.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_PRIMALDUALSYSTEMSTEP_H
11#define ROL_PRIMALDUALSYSTEMSTEP_H
12
15#include "ROL_SchurComplememt.hpp"
16
28namespace ROL {
29
30template<class Real>
31class PrimalDualSystemStep : public Step<Real> {
32
33 typedef Vector<Real> V;
40
45
46
47private:
48
49 // Block indices
50 static const size_type OPT = 0;
51 static const size_type EQUAL = 1;
52 static const size_type LOWER = 2;
53 static const size_type UPPER = 3;
54
55 // Super block indices
56 static const size_type OPTMULT = 0; // Optimization and equality multiplier components
57 static const size_type BNDMULT = 1; // Bound multiplier components
58
59 ROL::Ptr<Secant<Real> > secant_;
60 ROL::Ptr<Krylov<Real> > krylov_;
61 ROL::Ptr<V> scratch1_; // scratch vector
62 ROL::Ptr<V> scratch_;
63
64 ROL::Ptr<OP11> A_;
65 ROL::Ptr<OP12> B_;
66 ROL::Ptr<OP21> C_;
67 ROL::Ptr<OP22> D_;
68
69 ROL::Ptr<SCHUR> schur_; // Allows partial decoupling of (x,lambda) and (zl,zu)
70 ROL::Ptr<OP> op_; // Solve fully coupled system
71
75
78
79
80
81 // Repartition (x,lambda,zl,zu) as (xlambda,z) = ((x,lambda),(zl,zu))
82 ROL::Ptr<PV> repartition( V &x ) {
83
84 PV &x_pv = dynamic_cast<PV&>(x);
85 ROL::Ptr<V> xlambda = CreatePartitionedVector(x_pv.get(OPT),x_pv.get(EQUAL));
86 ROL::Ptr<V> z = CreatePartitionedVector(x_pv.get(LOWER),x_pv.get(UPPER));
87
88 ROL::Ptr<V> temp[] = {xlambda,z};
89
90 return ROL::makePtr<PV( std::vector<ROL::Ptr<V> >>(temp,temp+2) );
91
92 }
93
94 // Repartition (x,lambda,zl,zu) as (xlambda,z) = ((x,lambda),(zl,zu))
95 ROL::Ptr<const PV> repartition( const V &x ) {
96 const PV &x_pv = dynamic_cast<const PV&>(x);
97 ROL::Ptr<const V> xlambda = CreatePartitionedVector(x_pv.get(OPT),x_pv.get(EQUAL));
98 ROL::Ptr<const V> z = CreatePartitionedVector(x_pv.get(LOWER),x_pv.get(UPPER));
99
100 ROL::Ptr<const V> temp[] = {xlambda,z};
101
102 return ROL::makePtr<PV( std::vector<ROL::Ptr<const V> >>(temp,temp+2) );
103
104 }
105
106public:
107
108 using Step<Real>::initialize;
109 using Step<Real>::compute;
110 using Step<Real>::update;
111
112
113 PrimalDualSystemStep( ROL::ParameterList &parlist,
114 const ROL::Ptr<Krylov<Real> > &krylov,
115 const ROL::Ptr<Secant<Real> > &secant,
116 ROL::Ptr<V> &scratch1 ) : Step<Real>(),
117 krylov_(krylov), secant_(secant), scratch1_(scratch1), schur_(ROL::nullPtr),
118 op_(ROL::nullPtr), useSchurComplement_(false) {
119
120 PL &iplist = parlist.sublist("Step").sublist("Primal Dual Interior Point");
121 PL &syslist = iplist.sublist("System Solver");
122
123 useSchurComplement_ = syslist.get("Use Schur Complement",false);
124
125 }
126
127 PrimalDualSystemStep( ROL::ParameterList &parlist,
128 ROL::Ptr<V> &scratch1_ ) : Step<Real>() {
129 PrimalDualSystemStep(parlist,ROL::nullPtr,ROL::nullPtr,scratch1);
130 }
131
132 void initialize( V &x, const V &g, V &res, const V &c,
133 OBJ &obj, CON &con, BND &bnd, AS &algo_state ) {
134
135 Step<Real>::initialize(x,g,res,c,obj,con,bnd,algo_state);
136
137
138
139 ;
140
141 ROL::Ptr<OBJ> pObj = ROL::makePtrFromRef(obj);
142 ROL::Ptr<CON> pCon = ROL::makePtrFromRef(con);
143 ROL::Ptr<BND> pBnd = ROL::makePtrFromRef(bnd);
144
145 ROL::Ptr<PV> x_pv = repartition(x);
146
147 ROL::Ptr<V> xlambda = x_pv->get(OPTMULT);
148 ROL::Ptr<V> z = x_pv->get(BNDMULT);
149
150 A_ = ROL::makePtr<OP11>( pObj, pCon, *xlambda, scratch1_ );
151 B_ = ROL::makePtr<OP12>( );
152 C_ = ROL::makePtr<OP21>( *z );
153 D_ = ROL::makePtr<OP22>( pBnd, *xlambda );
154
155 if( useSchurComplement_ ) {
156 schur_ = ROL::makePtr<SCHUR>(A_,B_,C_,D_,scratch1_);
157 }
158 else {
160 }
161 }
162
163 void compute( V &s, const V &x, const V &res, OBJ &obj, CON &con,
164 BND &bnd, AS &algo_state ) {
165
166 ROL::Ptr<StepState<Real> > step_state = Step<Real>::getState();
167
168
169 if( useSchurComplement_ ) {
170
171 ROL::Ptr<const PV> x_pv = repartition(x);
172 ROL::Ptr<const PV> res_pv = repartition(res);
173 ROL::Ptr<PV> s_pv = repartition(s);
174
175
176 // Decouple (x,lambda) from (zl,zu) so that s <- L
177
178 ROL::Ptr<V> sxl = s_pv->get(OPTMULT);
179 ROL::Ptr<V> sz = s_pv->get(BNDMULT);
180
181
182
183 }
184 else {
185
186 }
187
188 }
189
190 void update( V &x, V &res, const V &s, OBJ &obj, CON &con,
191 BND &bnd, AS &algo_state ) {
192
193 ROL::Ptr<StepState<Real> > step_state = Step<Real>::getState();
194
195
196 }
197
198
199};
200
201} // namespace ROL
202
203#endif // ROL_PRIMALDUALSYSTEMSTEP_H
typename PV< Real >::size_type size_type
Provides the interface to apply a 2x2 block operator to a partitioned vector.
Provides the interface to apply upper and lower bound constraints.
Defines the general constraint operator interface.
Provides definitions for Krylov solvers.
Provides the interface to evaluate objective functions.
Defines the linear algebra of vector space on a generic partitioned vector.
ROL::Ptr< const Vector< Real > > get(size_type i) const
Provides the interface to compute approximate solutions to 2x2 block systems arising from primal-dual...
void initialize(V &x, const V &g, V &res, const V &c, OBJ &obj, CON &con, BND &bnd, AS &algo_state)
Initialize step with equality constraint.
PrimalDualInteriorPointBlock21 OP21
void update(V &x, V &res, const V &s, OBJ &obj, CON &con, BND &bnd, AS &algo_state)
Update step, if successful (equality constraints).
PrimalDualInteriorPointBlock22 OP22
int flagKrylov_
Termination flag for Krylov method (used for inexact Newton)
ROL::Ptr< const PV > repartition(const V &x)
ROL::Ptr< Krylov< Real > > krylov_
PrimalDualInteriorPointBlock11 OP11
ROL::Ptr< Secant< Real > > secant_
PrimalDualSystemStep(ROL::ParameterList &parlist, const ROL::Ptr< Krylov< Real > > &krylov, const ROL::Ptr< Secant< Real > > &secant, ROL::Ptr< V > &scratch1)
void compute(V &s, const V &x, const V &res, OBJ &obj, CON &con, BND &bnd, AS &algo_state)
Compute step (equality constraints).
int iterKrylov_
Number of Krylov iterations (used for inexact Newton)
PrimalDualSystemStep(ROL::ParameterList &parlist, ROL::Ptr< V > &scratch1_)
PrimalDualInteriorPointBlock12 OP12
Given a 2x2 block operator, perform the Schur reduction and return the decoupled system components.
Provides interface for and implements limited-memory secant operators.
Provides the interface to compute optimization steps.
Definition ROL_Step.hpp:34
virtual void initialize(Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &con, AlgorithmState< Real > &algo_state)
Initialize step with bound constraint.
Definition ROL_Step.hpp:54
ROL::Ptr< StepState< Real > > getState(void)
Definition ROL_Step.hpp:39
Defines the linear algebra or vector space interface.
ROL::Ptr< Vector< Real > > CreatePartitionedVector(const ROL::Ptr< Vector< Real > > &a)
State for algorithm class. Will be used for restarts.