ROL
poisson-control/example_01.cpp
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 
49 #define USE_HESSVEC 1
50 
51 #include "ROL_PoissonControl.hpp"
53 #include "ROL_TrustRegionStep.hpp"
54 #include "ROL_Algorithm.hpp"
55 #include "ROL_Types.hpp"
56 #include "Teuchos_oblackholestream.hpp"
57 #include "Teuchos_GlobalMPISession.hpp"
58 #include "Teuchos_XMLParameterListHelpers.hpp"
59 
60 #include <iostream>
61 #include <algorithm>
62 
63 #include "ROL_BoundConstraint.hpp"
64 
65 template<class Real>
67 private:
68  int dim_;
69  std::vector<Real> x_lo_;
70  std::vector<Real> x_up_;
71  Real min_diff_;
72 public:
73  BoundConstraint_PoissonControl(std::vector<Real> &lo, std::vector<Real> &up) {
74  dim_ = lo.size();
75  x_lo_.clear();
76  x_lo_.assign(lo.begin(),lo.end());
77  x_up_.clear();
78  x_up_.assign(up.begin(),up.end());
79  for ( unsigned i = 0; i < (unsigned)dim_; i++ ) {
80  if ( i == 0 ) {
81  min_diff_ = x_up_[i]-x_lo_[i];
82  }
83  else {
84  min_diff_ = std::min(min_diff_,x_up_[i]-x_lo_[i]);
85  }
86  }
87  min_diff_ *= 0.5;
88  }
89  bool isFeasible( const ROL::Vector<Real> &x ) {
90  Teuchos::RCP<const std::vector<Real> > ex =
91  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(x))).getVector();
92  bool val = true;
93  int cnt = 1;
94  for ( int i = 0; i < this->dim_; i++ ) {
95  if ( (*ex)[i] >= this->x_lo_[i] && (*ex)[i] <= this->x_up_[i] ) { cnt *= 1; }
96  else { cnt *= 0; }
97  }
98  if ( cnt == 0 ) { val = false; }
99  return val;
100  }
102  Teuchos::RCP<std::vector<Real> > ex =
103  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(x)).getVector());
104  for ( int i = 0; i < this->dim_; i++ ) {
105  (*ex)[i] = std::max(this->x_lo_[i],std::min(this->x_up_[i],(*ex)[i]));
106  }
107  }
108  void pruneLowerActive(ROL::Vector<Real> &v, const ROL::Vector<Real> &x, Real eps) {
109  Teuchos::RCP<const std::vector<Real> > ex =
110  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(x))).getVector();
111  Teuchos::RCP<std::vector<Real> > ev =
112  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(v)).getVector());
113  Real epsn = std::min(eps,this->min_diff_);
114  for ( int i = 0; i < this->dim_; i++ ) {
115  if ( ((*ex)[i] <= this->x_lo_[i]+epsn) ) {
116  (*ev)[i] = 0.0;
117  }
118  }
119  }
120  void pruneUpperActive(ROL::Vector<Real> &v, const ROL::Vector<Real> &x, Real eps) {
121  Teuchos::RCP<const std::vector<Real> > ex =
122  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(x))).getVector();
123  Teuchos::RCP<std::vector<Real> > ev =
124  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(v)).getVector());
125  Real epsn = std::min(eps,this->min_diff_);
126  for ( int i = 0; i < this->dim_; i++ ) {
127  if ( ((*ex)[i] >= this->x_up_[i]-epsn) ) {
128  (*ev)[i] = 0.0;
129  }
130  }
131  }
133  Teuchos::RCP<const std::vector<Real> > ex =
134  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(x))).getVector();
135  Teuchos::RCP<const std::vector<Real> > eg =
136  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(g))).getVector();
137  Teuchos::RCP<std::vector<Real> > ev =
138  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(v)).getVector());
139  Real epsn = std::min(eps,this->min_diff_);
140  for ( int i = 0; i < this->dim_; i++ ) {
141  if ( ((*ex)[i] <= this->x_lo_[i]+epsn && (*eg)[i] > 0.0) ){
142  (*ev)[i] = 0.0;
143  }
144  }
145  }
147  Teuchos::RCP<const std::vector<Real> > ex =
148  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(x))).getVector();
149  Teuchos::RCP<const std::vector<Real> > eg =
150  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(g))).getVector();
151  Teuchos::RCP<std::vector<Real> > ev =
152  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(v)).getVector());
153  Real epsn = std::min(eps,this->min_diff_);
154  for ( int i = 0; i < this->dim_; i++ ) {
155  if ( ((*ex)[i] >= this->x_up_[i]-epsn && (*eg)[i] < 0.0) ) {
156  (*ev)[i] = 0.0;
157  }
158  }
159  }
161  Teuchos::RCP<std::vector<Real> > us = Teuchos::rcp( new std::vector<Real>(this->dim_,0.0) );
162  us->assign(this->x_up_.begin(),this->x_up_.end());
163  Teuchos::RCP<ROL::Vector<Real> > up = Teuchos::rcp( new ROL::StdVector<Real>(us) );
164  u.set(*up);
165  }
167  Teuchos::RCP<std::vector<Real> > ls = Teuchos::rcp( new std::vector<Real>(this->dim_,0.0) );
168  ls->assign(this->x_lo_.begin(),this->x_lo_.end());
169  Teuchos::RCP<ROL::Vector<Real> > lp = Teuchos::rcp( new ROL::StdVector<Real>(ls) );
170  l.set(*lp);
171  }
172 // void pruneActive(ROL::Vector<Real> &v, const ROL::Vector<Real> &x, Real eps) {
173 // Teuchos::RCP<const std::vector<Real> > ex =
174 // (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(x))).getVector();
175 // Teuchos::RCP<std::vector<Real> > ev =
176 // Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(v)).getVector());
177 // Real epsn = std::min(eps,this->min_diff_);
178 // for ( int i = 0; i < this->dim_; i++ ) {
179 // if ( ((*ex)[i] <= this->x_lo_[i]+epsn) ||
180 // ((*ex)[i] >= this->x_up_[i]-epsn) ) {
181 // (*ev)[i] = 0.0;
182 // }
183 // }
184 // }
185 // void pruneActive(ROL::Vector<Real> &v, const ROL::Vector<Real> &g, const ROL::Vector<Real> &x, Real eps) {
186 // Teuchos::RCP<const std::vector<Real> > ex =
187 // (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(x))).getVector();
188 // Teuchos::RCP<const std::vector<Real> > eg =
189 // (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(g))).getVector();
190 // Teuchos::RCP<std::vector<Real> > ev =
191 // Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(v)).getVector());
192 // Real epsn = std::min(eps,this->min_diff_);
193 // for ( int i = 0; i < this->dim_; i++ ) {
194 // if ( ((*ex)[i] <= this->x_lo_[i]+epsn && (*eg)[i] > 0.0) ||
195 // ((*ex)[i] >= this->x_up_[i]-epsn && (*eg)[i] < 0.0) ) {
196 // (*ev)[i] = 0.0;
197 // }
198 // }
199 // }
200 };
201 
202 template <class Real>
203 class StatusTest_PDAS : public ROL::StatusTest<Real> {
204 private:
205 
206  Real gtol_;
207  Real stol_;
209 
210 public:
211 
212  virtual ~StatusTest_PDAS() {}
213 
214  StatusTest_PDAS( Real gtol = 1.e-6, Real stol = 1.e-12, int max_iter = 100 ) :
215  gtol_(gtol), stol_(stol), max_iter_(max_iter) {}
216 
219  virtual bool check( ROL::AlgorithmState<Real> &state ) {
220  if ( (state.gnorm > this->gtol_) &&
221  (state.snorm > this->stol_) &&
222  (state.iter < this->max_iter_) ) {
223  return true;
224  }
225  else {
226  if ( state.iter < 2 ) {
227  return true;
228  }
229  return false;
230  }
231  }
232 
233 };
234 
235 typedef double RealT;
236 
237 int main(int argc, char *argv[]) {
238 
239  Teuchos::GlobalMPISession mpiSession(&argc, &argv);
240 
241  // This little trick lets us print to std::cout only if a (dummy) command-line argument is provided.
242  int iprint = argc - 1;
243  Teuchos::RCP<std::ostream> outStream;
244  Teuchos::oblackholestream bhs; // outputs nothing
245  if (iprint > 0)
246  outStream = Teuchos::rcp(&std::cout, false);
247  else
248  outStream = Teuchos::rcp(&bhs, false);
249 
250  int errorFlag = 0;
251 
252  // *** Example body.
253 
254  try {
255  int dim = 256; // Set problem dimension.
256  RealT alpha = 1.e-4;
258  std::vector<RealT> lo(dim);
259  std::vector<RealT> up(dim);
260  for ( unsigned i = 0; i < (unsigned)dim; i++ ) {
261  if ( i < (unsigned)dim/3 || i > 2*(unsigned)dim/3 ) {
262  lo[i] = 0.0;
263  up[i] = 0.25;
264  }
265  else {
266  lo[i] = 0.75;
267  up[i] = 1.0;
268  }
269  }
271 
272  Teuchos::ParameterList parlist;
273  // Krylov parameters.
274  parlist.set("Absolute Krylov Tolerance", 1.e-4);
275  parlist.set("Relative Krylov Tolerance", 1.e-2);
276  parlist.set("Maximum Number of Krylov Iterations", 50);
277  // PDAS parameters.
278  parlist.set("PDAS Relative Step Tolerance", 1.e-8);
279  parlist.set("PDAS Relative Gradient Tolerance", 1.e-6);
280  parlist.set("PDAS Maximum Number of Iterations", 1);
281  parlist.set("PDAS Dual Scaling", alpha);
282  // Define step.
283  ROL::PrimalDualActiveSetStep<RealT> step_pdas(parlist);
284 
285  // Define status test.
286  RealT gtol = 1e-12; // norm of gradient tolerance
287  RealT stol = 1e-14; // norm of step tolerance
288  int maxit = 100; // maximum number of iterations
289  StatusTest_PDAS<RealT> status_pdas(gtol, stol, maxit);
290 
291  // Define algorithm.
292  ROL::DefaultAlgorithm<RealT> algo_pdas(step_pdas,status_pdas,false);
293 
294  // Iteration vector.
295  Teuchos::RCP<std::vector<RealT> > x_rcp = Teuchos::rcp( new std::vector<RealT> (dim, 0.0) );
296  // Set initial guess.
297  for (int i=0; i<dim; i++) {
298  (*x_rcp)[i] = 0.0;
299  }
300  ROL::StdVector<RealT> x(x_rcp);
301 
302  // Run algorithm.
303  algo_pdas.run(x, obj, icon, true);
304 
305  std::ofstream file;
306  file.open("control_PDAS.txt");
307  for ( unsigned i = 0; i < (unsigned)dim; i++ ) {
308  file << (*x_rcp)[i] << "\n";
309  }
310  file.close();
311 
312  // Use Projected Methods
313  std::string filename = "input.xml";
314  Teuchos::RCP<Teuchos::ParameterList> parlist_tr = Teuchos::rcp( new Teuchos::ParameterList() );
315  Teuchos::updateParametersFromXmlFile( filename, Teuchos::Ptr<Teuchos::ParameterList>(&*parlist_tr) );
316  ROL::TrustRegionStep<RealT> step_tr(*parlist_tr);
317  ROL::StatusTest<RealT> status_tr(gtol,stol,maxit);
318  ROL::DefaultAlgorithm<RealT> algo_tr(step_tr,status_tr,false);
319  // Iteration vector.
320  Teuchos::RCP<std::vector<RealT> > y_rcp = Teuchos::rcp( new std::vector<RealT> (dim, 0.0) );
321  // Set initial guess.
322  for (int i=0; i<dim; i++) {
323  (*y_rcp)[i] = 0.0;
324  }
325  ROL::StdVector<RealT> y(y_rcp);
326  // Run Algorithm
327  algo_tr.run(y,obj,icon,true);
328 
329  std::ofstream file_tr;
330  file_tr.open("control_TR.txt");
331  for ( unsigned i = 0; i < (unsigned)dim; i++ ) {
332  file_tr << (*y_rcp)[i] << "\n";
333  }
334  file_tr.close();
335 
336  Teuchos::RCP<ROL::Vector<RealT> > error = x.clone();
337  error->set(x);
338  error->axpy(-1.0,y);
339  std::cout << "\nError between PDAS solution and TR solution is " << error->norm() << "\n";
340 
341  }
342  catch (std::logic_error err) {
343  *outStream << err.what() << "\n";
344  errorFlag = -1000;
345  }; // end try
346 
347  if (errorFlag != 0)
348  std::cout << "End Result: TEST FAILED\n";
349  else
350  std::cout << "End Result: TEST PASSED\n";
351 
352  return 0;
353 
354 }
355 
Implements the computation of optimization steps with the Newton primal-dual active set method...
int main(int argc, char *argv[])
StatusTest_PDAS(Real gtol=1.e-6, Real stol=1.e-12, int max_iter=100)
void project(ROL::Vector< Real > &x)
Project optimization variables onto the bounds.
BoundConstraint_PoissonControl(std::vector< Real > &lo, std::vector< Real > &up)
Contains definitions of custom data types in ROL.
Contains definitions for Poisson optimal control.
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:72
Teuchos::RCP< Vector< Real > > clone() const
Clone to make a new (uninitialized) vector.
State for algorithm class. Will be used for restarts.
Definition: ROL_Types.hpp:76
void setVectorToLowerBound(ROL::Vector< Real > &l)
Set the input vector to the lower bound.
Provides the std::vector implementation of the ROL::Vector interface.
virtual std::vector< std::string > run(Vector< Real > &x, Objective< Real > &obj, bool print=false, std::ostream &outStream=std::cout)
Run algorithm on unconstrained problems (Type-U). This is the primary Type-U interface.
Poisson distributed control.
Provides an interface to check status of optimization algorithms.
void pruneUpperActive(ROL::Vector< Real > &v, const ROL::Vector< Real > &x, Real eps)
Set variables to zero if they correspond to the upper -active set.
bool isFeasible(const ROL::Vector< Real > &x)
Check if the vector, v, is feasible.
Provides the interface to apply upper and lower bound constraints.
void pruneLowerActive(ROL::Vector< Real > &v, const ROL::Vector< Real > &g, const ROL::Vector< Real > &x, Real eps)
Set variables to zero if they correspond to the lower -binding set.
void pruneLowerActive(ROL::Vector< Real > &v, const ROL::Vector< Real > &x, Real eps)
Set variables to zero if they correspond to the lower -active set.
double RealT
virtual void set(const Vector &x)
Set where .
Definition: ROL_Vector.hpp:194
virtual bool check(ROL::AlgorithmState< Real > &state)
Check algorithm status.
void setVectorToUpperBound(ROL::Vector< Real > &u)
Set the input vector to the upper bound.
void pruneUpperActive(ROL::Vector< Real > &v, const ROL::Vector< Real > &g, const ROL::Vector< Real > &x, Real eps)
Set variables to zero if they correspond to the upper -binding set.
Provides the interface to compute optimization steps with trust regions.