DOpE
reducednewtonalgorithm.h
Go to the documentation of this file.
1 
24 #ifndef _REDUCEDNEWTON__ALGORITHM_H_
25 #define _REDUCEDNEWTON__ALGORITHM_H_
26 
27 #include "reducedalgorithm.h"
28 #include "parameterreader.h"
29 
30 #include <iostream>
31 #include <assert.h>
32 #include <iomanip>
33 namespace DOpE
34 {
49  template <typename PROBLEM, typename VECTOR>
50  class ReducedNewtonAlgorithm : public ReducedAlgorithm<PROBLEM, VECTOR>
51  {
52  public:
66  ReducedNewtonAlgorithm(PROBLEM* OP,
68  ParameterReader &param_reader,
69  DOpEExceptionHandler<VECTOR>* Except=NULL,
70  DOpEOutputHandler<VECTOR>* Output=NULL,
71  int base_priority=0);
73 
79  static void declare_params(ParameterReader &param_reader);
80 
90  virtual int Solve(ControlVector<VECTOR>& q,double global_tol=-1.);
97  double NewtonResidual(const ControlVector<VECTOR>& q);
98 
99  protected:
118  virtual int SolveReducedLinearSystem(const ControlVector<VECTOR>& q,
119  const ControlVector<VECTOR>& gradient,
120  const ControlVector<VECTOR>& gradient_transposed,
122 
136  virtual int ReducedNewtonLineSearch(const ControlVector<VECTOR>& dq,
137  const ControlVector<VECTOR>& gradient,
138  double& cost,
143  virtual double Residual(const ControlVector<VECTOR>& gradient,
144  const ControlVector<VECTOR>& gradient_transposed)
145  {return gradient*gradient_transposed;}
146  private:
147  unsigned int _nonlinear_maxiter, _linear_maxiter, _line_maxiter;
148  double _nonlinear_tol, _nonlinear_global_tol, _linear_tol, _linear_global_tol, _linesearch_rho, _linesearch_c;
149  bool _compute_functionals_in_every_step;
150  std::string _postindex;
151  };
152 
153  /***************************************************************************************/
154  /****************************************IMPLEMENTATION*********************************/
155  /***************************************************************************************/
156  using namespace dealii;
157 
158  /******************************************************/
159 
160 template <typename PROBLEM, typename VECTOR>
162  {
163  param_reader.SetSubsection("reducednewtonalgorithm parameters");
164  param_reader.declare_entry("nonlinear_maxiter", "10",Patterns::Integer(0));
165  param_reader.declare_entry("nonlinear_tol", "1.e-7",Patterns::Double(0));
166  param_reader.declare_entry("nonlinear_global_tol", "1.e-11",Patterns::Double(0));
167 
168  param_reader.declare_entry("linear_maxiter", "40",Patterns::Integer(0));
169  param_reader.declare_entry("linear_tol", "1.e-10",Patterns::Double(0));
170  param_reader.declare_entry("linear_global_tol", "1.e-12",Patterns::Double(0));
171 
172  param_reader.declare_entry("line_maxiter", "4",Patterns::Integer(0));
173  param_reader.declare_entry("linesearch_rho", "0.9",Patterns::Double(0));
174  param_reader.declare_entry("linesearch_c", "0.1",Patterns::Double(0));
175 
176  param_reader.declare_entry("compute_functionals_in_every_step", "false",Patterns::Bool());
177 
179  }
180 /******************************************************/
181 
182 template <typename PROBLEM, typename VECTOR>
185  ParameterReader &param_reader,
188  int base_priority)
189  : ReducedAlgorithm<PROBLEM, VECTOR>(OP,S,param_reader,Except,Output,base_priority)
190  {
191 
192  param_reader.SetSubsection("reducednewtonalgorithm parameters");
193  _nonlinear_maxiter = param_reader.get_integer ("nonlinear_maxiter");
194  _nonlinear_tol = param_reader.get_double ("nonlinear_tol");
195  _nonlinear_global_tol = param_reader.get_double ("nonlinear_global_tol");
196 
197  _linear_maxiter = param_reader.get_integer ("linear_maxiter");
198  _linear_tol = param_reader.get_double ("linear_tol");
199  _linear_global_tol = param_reader.get_double ("linear_global_tol");
200 
201  _line_maxiter = param_reader.get_integer ("line_maxiter");
202  _linesearch_rho = param_reader.get_double ("linesearch_rho");
203  _linesearch_c = param_reader.get_double ("linesearch_c");
204 
205  _compute_functionals_in_every_step = param_reader.get_bool ("compute_functionals_in_every_step");
206 
207  _postindex = "_"+this->GetProblem()->GetName();
208  }
209 
210 /******************************************************/
211 
212 template <typename PROBLEM, typename VECTOR>
214  {
215 
216  }
217 
218 /******************************************************/
219 
220 template <typename PROBLEM, typename VECTOR>
222 {
223  //Solve j'(q) = 0
224  ControlVector<VECTOR> gradient(q), gradient_transposed(q);
225 
226  try
227  {
228  this->GetReducedProblem()->ComputeReducedCostFunctional(q);
229  }
230  catch(DOpEException& e)
231  {
232  this->GetExceptionHandler()->HandleCriticalException(e,"ReducedNewtonAlgorithm::NewtonResidual");
233  }
234 
235  try
236  {
237  this->GetReducedProblem()->ComputeReducedGradient(q,gradient,gradient_transposed);
238  }
239  catch(DOpEException& e)
240  {
241  this->GetExceptionHandler()->HandleCriticalException(e,"ReducedNewtonAlgorithm::NewtonResidual");
242  }
243 
244  return sqrt(Residual(gradient,gradient_transposed));
245 }
246 
247 /******************************************************/
248 
249 template <typename PROBLEM, typename VECTOR>
251 {
252 
253  q.ReInit();
254  //Solve j'(q) = 0
255  ControlVector<VECTOR> dq(q), gradient(q), gradient_transposed(q);
256 
257  unsigned int iter=0;
258  double cost=0.;
259  std::stringstream out;
260  this->GetOutputHandler()->InitNewtonOut(out);
261 
262  out << "**************************************************\n";
263  out << "* Starting Reduced Newton Algorithm *\n";
264  out << "* Solving : "<<this->GetProblem()->GetName()<<"\t*\n";
265  out << "* CDoFs : ";
266  q.PrintInfos(out);
267  out << "* SDoFs : ";
268  this->GetReducedProblem()->StateSizeInfo(out);
269  out << "**************************************************";
270  this->GetOutputHandler()->Write(out,1+this->GetBasePriority(),1,1);
271 
272  this->GetOutputHandler()->SetIterationNumber(iter,"OptNewton"+_postindex);
273 
274  this->GetOutputHandler()->Write(q,"Control"+_postindex,"control");
275 
276  try
277  {
278  cost = this->GetReducedProblem()->ComputeReducedCostFunctional(q);
279  }
280  catch(DOpEException& e)
281  {
282  this->GetExceptionHandler()->HandleCriticalException(e,"ReducedNewtonAlgorithm::Solve");
283  }
284 
285  out<< "CostFunctional: " << cost;
286  this->GetOutputHandler()->Write(out,2+this->GetBasePriority());
287 
288  if (_compute_functionals_in_every_step == true)
289  {
290  try
291  {
292  this->GetReducedProblem()->ComputeReducedFunctionals(q);
293  }
294  catch(DOpEException& e)
295  {
296  this->GetExceptionHandler()->HandleCriticalException(e);
297  }
298  }
299 
300  try
301  {
302  this->GetReducedProblem()->ComputeReducedGradient(q,gradient,gradient_transposed);
303  }
304  catch(DOpEException& e)
305  {
306  this->GetExceptionHandler()->HandleCriticalException(e,"ReducedNewtonAlgorithm::Solve");
307  }
308 
309  double res = Residual(gradient,gradient_transposed);//gradient*gradient_transposed;
310  double firstres = res;
311 
312  assert(res >= 0);
313 
314  this->GetOutputHandler()->Write(gradient,"NewtonResidual"+_postindex,"control");
315  out<< "\t Newton step: " <<iter<<"\t Residual (abs.): "<<sqrt(res)<<"\n";
316  out<< "\t Newton step: " <<iter<<"\t Residual (rel.): "<<std::scientific<<sqrt(res)/sqrt(res)<<"\n";
317  this->GetOutputHandler()->Write(out,3+this->GetBasePriority());
318  int liniter = 0;
319  int lineiter =0;
320  unsigned int miniter = 0;
321  if(global_tol > 0.)
322  miniter = 1;
323 
324  global_tol = std::max(_nonlinear_global_tol,global_tol);
325  while(( (res >= global_tol*global_tol) && (res >= _nonlinear_tol*_nonlinear_tol*firstres) ) || iter < miniter )
326  {
327  iter++;
328  this->GetOutputHandler()->SetIterationNumber(iter,"OptNewton"+_postindex);
329 
330  if(iter > _nonlinear_maxiter)
331  {
332  throw DOpEIterationException("Iteration count exceeded bounds!","ReducedNewtonAlgorithm::Solve");
333  }
334 
335  //Compute a search direction
336  try
337  {
338  liniter = SolveReducedLinearSystem(q,gradient,gradient_transposed, dq);
339  }
340  catch(DOpEIterationException& e)
341  {
342  //Seems uncritical too many linear solves, it'll probably work
343  //So only write a warning, and continue.
344  this->GetExceptionHandler()->HandleException(e,"ReducedNewtonAlgorithm::Solve");
345  liniter = -1;
346  }
348  {
349  this->GetExceptionHandler()->HandleException(e,"ReducedNewtonAlgorithm::Solve");
350  lineiter = -2;
351  }
352  catch(DOpEException& e)
353  {
354  this->GetExceptionHandler()->HandleCriticalException(e,"ReducedNewtonAlgorithm::Solve");
355  }
356  //Linesearch
357  try
358  {
359  lineiter = ReducedNewtonLineSearch(dq,gradient,cost,q);
360  }
361  catch(DOpEIterationException& e)
362  {
363  //Seems uncritical too many line search steps, it'll probably work
364  //So only write a warning, and continue.
365  this->GetExceptionHandler()->HandleException(e,"ReducedNewtonAlgorithm::Solve");
366  lineiter = -1;
367  }
368  //catch(DOpEException& e)
369  //{
370  // this->GetExceptionHandler()->HandleCriticalException(e);
371  //}
372 
373  out<< "CostFunctional: " << cost;
374  this->GetOutputHandler()->Write(out,3+this->GetBasePriority());
375 
376  if (_compute_functionals_in_every_step == true)
377  {
378  try
379  {
380  this->GetReducedProblem()->ComputeReducedFunctionals(q);
381  }
382  catch(DOpEException& e)
383  {
384  this->GetExceptionHandler()->HandleCriticalException(e);
385  }
386  }
387 
388 
389  //Prepare the next Iteration
390  try
391  {
392  this->GetReducedProblem()->ComputeReducedGradient(q,gradient,gradient_transposed);
393  }
394  catch(DOpEException& e)
395  {
396  this->GetExceptionHandler()->HandleCriticalException(e,"ReducedNewtonAlgorithm::Solve");
397  }
398 
399  this->GetOutputHandler()->Write(q,"Control"+_postindex,"control");
400  this->GetOutputHandler()->Write(gradient,"NewtonResidual"+_postindex,"control");
401 
402  res = Residual(gradient,gradient_transposed);//gradient*gradient_transposed;
403 
404  out<<"\t Newton step: " <<iter<<"\t Residual (rel.): "<<this->GetOutputHandler()->ZeroTolerance(sqrt(res)/sqrt(firstres),1.0)<< "\t LinearIters ["<<liniter<<"]\t LineSearch {"<<lineiter<<"} ";
405  this->GetOutputHandler()->Write(out,3+this->GetBasePriority());
406  }
407 
408  //We are done write total evaluation
409  out<< "CostFunctional: " << cost;
410  this->GetOutputHandler()->Write(out,2+this->GetBasePriority());
411  try
412  {
413  this->GetReducedProblem()->ComputeReducedFunctionals(q);
414  }
415  catch(DOpEException& e)
416  {
417  this->GetExceptionHandler()->HandleCriticalException(e,"ReducedNewtonAlgorithm::Solve");
418  }
419 
420  out << "**************************************************\n";
421  out << "* Stopping Reduced Newton Algorithm *\n";
422  out << "* after "<<std::setw(6)<<iter<<" Iterations *\n";
423  out.precision(4);
424  out << "* with rel. Residual "<<std::scientific << std::setw(11) << this->GetOutputHandler()->ZeroTolerance(sqrt(res)/sqrt(firstres),1.0)<<" *\n";
425  out.precision(10);
426  out << "**************************************************";
427  this->GetOutputHandler()->Write(out,1+this->GetBasePriority(),1,1);
428  return iter;
429 }
430 
431 /******************************************************/
432 
433 template <typename PROBLEM, typename VECTOR>
435  const ControlVector<VECTOR>& gradient,
436  const ControlVector<VECTOR>& gradient_transposed,
438 {
439  std::stringstream out;
440  dq = 0.;
441  ControlVector<VECTOR> r(q), r_transposed(q), d(q), Hd(q), Hd_transposed(q);
442 
443  r = gradient;
444  r_transposed = gradient_transposed;
445  d = gradient_transposed;
446 
447  double res = Residual(r,r_transposed);//r*r_transposed;
448  double firstres = res;
449 
450  assert(res >= 0.);
451 
452  out << "Starting Reduced Linear Solver with Residual: "<<sqrt(res);
453  this->GetOutputHandler()->Write(out,4+this->GetBasePriority());
454 
455  unsigned int iter = 0;
456  double cgalpha, cgbeta, oldres;
457 
458  this->GetOutputHandler()->SetIterationNumber(iter,"OptNewtonCg"+_postindex);
459 
460  //while(res>=_linear_tol*_linear_tol*firstres && res>=_linear_global_tol*_linear_global_tol)
461  //using Algorithm 6.1 from Nocedal Wright
462  while(res>= std::min(0.25,sqrt(firstres))*firstres && res>=_linear_global_tol*_linear_global_tol)
463  {
464  iter++;
465  this->GetOutputHandler()->SetIterationNumber(iter,"OptNewtonCg"+_postindex);
466  if(iter > _linear_maxiter)
467  {
468  throw DOpEIterationException("Iteration count exceeded bounds!","ReducedNewtonAlgorithm::SolveReducedLinearSystem");
469  }
470 
471  try
472  {
473  this->GetReducedProblem()->ComputeReducedHessianVector(q,d,Hd,Hd_transposed);
474  }
475  catch(DOpEException& e)
476  {
477  this->GetExceptionHandler()->HandleCriticalException(e,"ReducedNewtonAlgorithm::SolveReducedLinearSystem");
478  }
479 
480  cgalpha = res / (Hd*d);
481 
482  if(cgalpha < 0)
483  {
484  if(iter==1)
485  {
486  dq.add(cgalpha,d);
487  }
488  throw DOpENegativeCurvatureException("Negative curvature detected!","ReducedNewtonAlgorithm::SolveReducedLinearSystem");
489  }
490 
491  dq.add(cgalpha,d);
492  r.add(cgalpha,Hd);
493  r_transposed.add(cgalpha,Hd_transposed);
494 
495  oldres = res;
496  res = Residual(r,r_transposed);//r*r_transposed;
497  if(res < 0.)
498  {
499  //something is broken, maybe don't use update formula and
500  //calculate res from scratch.
501  try
502  {
503  this->GetReducedProblem()->ComputeReducedHessianVector(q,dq,Hd,Hd_transposed);
504  }
505  catch(DOpEException& e)
506  {
507  this->GetExceptionHandler()->HandleCriticalException(e);
508  }
509  r = gradient;
510  r_transposed = gradient_transposed;
511  r.add(1.,Hd);
512  r_transposed.add(1.,Hd_transposed);
513  res = Residual(r,r_transposed);
514  }
515  assert(res >= 0.);
516  out<<"\t Cg step: " <<iter<<"\t Residual: "<<sqrt(res);
517  this->GetOutputHandler()->Write(out,4+this->GetBasePriority());
518 
519  cgbeta = res / oldres; //Fletcher-Reeves
520  d*= cgbeta;
521  d.equ(-1,r_transposed);
522  }
523  return iter;
524 }
525 
526 
527 /******************************************************/
528 
529 template <typename PROBLEM, typename VECTOR>
531  const ControlVector<VECTOR>& gradient,
532  double& cost,
534 {
535  double rho = _linesearch_rho;
536  double c = _linesearch_c;
537 
538  double costnew = 0.;
539  bool force_linesearch = false;
540 
541  q+=dq;
542  try
543  {
544  costnew = this->GetReducedProblem()->ComputeReducedCostFunctional(q);
545  }
546  catch(DOpEException& e)
547  {
548 // this->GetExceptionHandler()->HandleException(e);
549  force_linesearch = true;
550  this->GetOutputHandler()->Write("Computing Cost Failed",4+this->GetBasePriority());
551  }
552 
553  double alpha=1;
554  unsigned int iter =0;
555 
556  double reduction = gradient*dq;
557  if(reduction > 0)
558  {
559  this->GetOutputHandler()->WriteError("Waring: computed direction doesn't seem to be a descend direction!");
560  reduction = 0;
561  }
562 
563  if(_line_maxiter > 0)
564  {
565  if(fabs(reduction) < 1.e-10*cost)
566  reduction = 0.;
567  if(std::isinf(costnew) || std::isnan(costnew) || (costnew >= cost + c*alpha*reduction) || force_linesearch)
568  {
569  this->GetOutputHandler()->Write("\t linesearch ",4+this->GetBasePriority());
570  while(std::isinf(costnew) || std::isnan(costnew) || (costnew >= cost + c*alpha*reduction) || force_linesearch)
571  {
572  iter++;
573  if(iter > _line_maxiter)
574  {
575  if(force_linesearch)
576  {
577  throw DOpEException("Iteration count exceeded bounds while unable to compute the CostFunctional!","ReducedNewtonAlgorithm::ReducedNewtonLineSearch");
578  }
579  else
580  {
581  cost = costnew;
582  throw DOpEIterationException("Iteration count exceeded bounds!","ReducedNewtonAlgorithm::ReducedNewtonLineSearch");
583  }
584  }
585  force_linesearch = false;
586  q.add(alpha*(rho-1.),dq);
587  alpha *= rho;
588 
589  try
590  {
591  costnew = this->GetReducedProblem()->ComputeReducedCostFunctional(q);
592  }
593  catch(DOpEException& e)
594  {
595  //this->GetExceptionHandler()->HandleException(e);
596  force_linesearch = true;
597  this->GetOutputHandler()->Write("Computing Cost Failed",4+this->GetBasePriority());
598  }
599  }
600  }
601  cost = costnew;
602  }
603 
604  return iter;
605 
606 }
607 
608 
609 }
610 #endif
Definition: dopeexception.h:76
void PrintInfos(std::stringstream &out)
Definition: controlvector.cc:542
ReducedNewtonAlgorithm(PROBLEM *OP, ReducedProblemInterface< PROBLEM, VECTOR > *S, ParameterReader &param_reader, DOpEExceptionHandler< VECTOR > *Except=NULL, DOpEOutputHandler< VECTOR > *Output=NULL, int base_priority=0)
Definition: reducednewtonalgorithm.h:183
double get_double(const std::string &entry_name)
Definition: parameterreader.h:115
Definition: dopeexception.h:90
bool get_bool(const std::string &entry_name)
Definition: parameterreader.h:148
Definition: parameterreader.h:36
Definition: optproblemcontainer.h:70
~ReducedNewtonAlgorithm()
Definition: reducednewtonalgorithm.h:213
virtual double Residual(const ControlVector< VECTOR > &gradient, const ControlVector< VECTOR > &gradient_transposed)
Definition: reducednewtonalgorithm.h:143
void declare_entry(const std::string &entry, const std::string &default_value, const Patterns::PatternBase &pattern=Patterns::Anything(), const std::string &documentation=std::string())
Definition: parameterreader.h:98
Definition: controlvector.h:48
Definition: reducednewtonalgorithm.h:50
PROBLEM * GetProblem()
Definition: reducedalgorithm.h:187
virtual int Solve(ControlVector< VECTOR > &q, double global_tol=-1.)
Definition: reducednewtonalgorithm.h:250
virtual int SolveReducedLinearSystem(const ControlVector< VECTOR > &q, const ControlVector< VECTOR > &gradient, const ControlVector< VECTOR > &gradient_transposed, ControlVector< VECTOR > &dq)
Definition: reducednewtonalgorithm.h:434
Definition: reducedalgorithm.h:54
void add(double s, const ControlVector &dq)
Definition: controlvector.cc:363
void ReInit()
Definition: controlvector.cc:60
static void declare_params(ParameterReader &param_reader)
Definition: reducedalgorithm.h:241
virtual int ReducedNewtonLineSearch(const ControlVector< VECTOR > &dq, const ControlVector< VECTOR > &gradient, double &cost, ControlVector< VECTOR > &q)
Definition: reducednewtonalgorithm.h:530
int get_integer(const std::string &entry_name)
Definition: parameterreader.h:126
static void declare_params(ParameterReader &param_reader)
Definition: reducednewtonalgorithm.h:161
void SetSubsection(const std::string subsection)
Definition: parameterreader.h:93
Definition: reducedprobleminterface.h:338
Definition: dopeexception.h:35
double NewtonResidual(const ControlVector< VECTOR > &q)
Definition: reducednewtonalgorithm.h:221
Definition: optproblemcontainer.h:72