queso-0.56.1
GslOptimizer.C
Go to the documentation of this file.
1 //-----------------------------------------------------------------------bl-
2 //--------------------------------------------------------------------------
3 //
4 // QUESO - a library to support the Quantification of Uncertainty
5 // for Estimation, Simulation and Optimization
6 //
7 // Copyright (C) 2008-2015 The PECOS Development Team
8 //
9 // This library is free software; you can redistribute it and/or
10 // modify it under the terms of the Version 2.1 GNU Lesser General
11 // Public License as published by the Free Software Foundation.
12 //
13 // This library is distributed in the hope that it will be useful,
14 // but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // Lesser General Public License for more details.
17 //
18 // You should have received a copy of the GNU Lesser General Public
19 // License along with this library; if not, write to the Free Software
20 // Foundation, Inc. 51 Franklin Street, Fifth Floor,
21 // Boston, MA 02110-1301 USA
22 //
23 //-----------------------------------------------------------------------el-
24 
25 #include <iostream>
26 
27 #include <queso/Defines.h>
28 #include <queso/GslVector.h>
29 #include <queso/VectorSpace.h>
30 #include <queso/ScalarFunction.h>
31 #include <queso/GslOptimizer.h>
32 #include <queso/OptimizerMonitor.h>
33 
34 #include <gsl/gsl_multimin.h>
35 #include <gsl/gsl_blas.h>
36 
37 namespace QUESO {
38 
39 // We need to extern "C" because gsl needs a pointer to a C function to
40 // minimize
41 extern "C" {
42  // This evaluate -log posterior
43  double c_evaluate(const gsl_vector * x, void * context) {
44 
45  GslOptimizer * optimizer = static_cast<GslOptimizer * >(context);
46 
47  GslVector state(
48  optimizer->objectiveFunction().domainSet().vectorSpace().zeroVector());
49 
50  // DM: Doing this copy sucks, but whatever. It'll do for now.
51  for (unsigned int i = 0; i < state.sizeLocal(); i++) {
52  state[i] = gsl_vector_get(x, i);
53  }
54 
55  // Bail early if GSL tries to evaluate outside of the domain
56  if (!optimizer->objectiveFunction().domainSet().contains(state)) {
57  return GSL_NAN;
58  }
59 
60  // Should cache derivative here so we don't a) segfault in the user's code
61  // and b) so we don't recompute stuff
62  double result = -optimizer->objectiveFunction().lnValue(state, NULL, NULL,
63  NULL, NULL);
64 
65  return result;
66  }
67 
68  // This evaluates the derivative of -log posterior
69  void c_evaluate_derivative(const gsl_vector * x, void * context,
70  gsl_vector * derivative) {
71  GslOptimizer * optimizer = static_cast<GslOptimizer * >(context);
72 
73  GslVector state(
74  optimizer->objectiveFunction().domainSet().vectorSpace().zeroVector());
75  GslVector deriv(
76  optimizer->objectiveFunction().domainSet().vectorSpace().zeroVector());
77 
78  // DM: Doing this copy sucks, but whatever. It'll do for now.
79  for (unsigned int i = 0; i < state.sizeLocal(); i++) {
80  state[i] = gsl_vector_get(x, i);
81 
82  // We fill with GSL_NAN and use it as a flag to check later that the user
83  // actually fills the derivative vector with stuff
84  deriv[i] = GSL_NAN;
85  }
86 
87  if (!optimizer->objectiveFunction().domainSet().contains(state)) {
88  // Fill derivative with error codes if the point is outside of the
89  // domain
90  for (unsigned int i = 0; i < deriv.sizeLocal(); i++) {
91  gsl_vector_set(derivative, i, GSL_NAN);
92  }
93  }
94  else {
95  // We should cache the return value here so we don't recompute stuff
96  double fx = -optimizer->objectiveFunction().lnValue(state, NULL, &deriv,
97  NULL, NULL);
98 
99  // Decide whether or not we need to do a finite difference based on
100  // whether the user actually filled deriv with values that are not
101  // GSL_NAN
102  //
103  // We're currently doing this check every time this function gets called.
104  // We could probably pull this logic out of here and put it somewhere
105  // where it only happens once
106  bool userComputedDerivative = true;
107  for (unsigned int i = 0; i < deriv.sizeLocal(); i++) {
108  // If the user missed out a derivative in any direction, fall back to
109  // a finite difference
110  if (gsl_isnan(deriv[i])) {
111  userComputedDerivative = false;
112  break;
113  }
114  }
115 
116  if (userComputedDerivative) {
117  for (unsigned int i = 0; i < deriv.sizeLocal(); i++) {
118  gsl_vector_set(derivative, i, -deriv[i]); // We need the minus sign
119  }
120  }
121  else {
122  // Finite difference step-size
123  double h = optimizer->getFiniteDifferenceStepSize();
124 
125  // User did not provide a derivative, so do a finite difference
126  for (unsigned int i = 0; i < deriv.sizeLocal(); i++) {
127  double tempState = state[i];
128  state[i] += h;
129 
130  // User didn't provide a derivative, so we don't bother passing in
131  // the derivative vector again
132  double fxph = -optimizer->objectiveFunction().lnValue(state, NULL,
133  NULL, NULL, NULL);
134 
135  // Reset the state back to what it was before
136  state[i] = tempState;
137 
138  // Make sure we didn't do anything dumb and tell gsl if we did
139  if (!gsl_isnan(fx) && !gsl_isnan(fxph)) {
140  gsl_vector_set(derivative, i, (fxph - fx) / h);
141  }
142  else {
143  gsl_vector_set(derivative, i, GSL_NAN);
144  }
145  }
146  }
147  }
148  }
149 
150  // This evaluates -log posterior and the derivative of -log posterior
151  void c_evaluate_with_derivative(const gsl_vector * x, void * context,
152  double * f, gsl_vector * derivative) {
153  // We don't need to call both of these
154  *f = c_evaluate(x, context);
155  c_evaluate_derivative(x, context, derivative);
156  }
157 } // End extern "C"
158 
160  const BaseScalarFunction<GslVector, GslMatrix> & objectiveFunction)
161  : BaseOptimizer(),
162  m_objectiveFunction(objectiveFunction),
163  m_initialPoint(new GslVector(objectiveFunction.domainSet().
164  vectorSpace().zeroVector())),
165  m_minimizer(new GslVector(this->m_objectiveFunction.domainSet().
166  vectorSpace().zeroVector())),
167  m_solver_type(BFGS2),
168  m_fstep_size(this->m_objectiveFunction.domainSet().vectorSpace().zeroVector()),
169  m_fdfstep_size(getFdfstepSize()),
170  m_line_tol(getLineTolerance())
171 {
172  // We initialize the minimizer to GSL_NAN just in case the optimization fails
173  m_minimizer->cwSet(GSL_NAN);
174 
175  // Set to documented default value.
177 
178  // Set solver type to the one set in the options object
180 }
181 
183  OptimizerOptions options,
184  const BaseScalarFunction<GslVector, GslMatrix> & objectiveFunction)
185  : BaseOptimizer(options),
186  m_objectiveFunction(objectiveFunction),
187  m_initialPoint(new GslVector(objectiveFunction.domainSet().
188  vectorSpace().zeroVector())),
189  m_minimizer(new GslVector(this->m_objectiveFunction.domainSet().
190  vectorSpace().zeroVector())),
191  m_solver_type(BFGS2),
192  m_fstep_size(this->m_objectiveFunction.domainSet().vectorSpace().zeroVector()),
193  m_fdfstep_size(getFdfstepSize()),
194  m_line_tol(getLineTolerance())
195 {
196  // We initialize the minimizer to GSL_NAN just in case the optimization fails
197  m_minimizer->cwSet(GSL_NAN);
198 
199  // Set to documented default value.
201 
202  // Set solver type to the one set in the options object
204 }
205 
207 {
208  delete this->m_initialPoint;
209 }
210 
212 
213  // First check that initial guess is reasonable
214  if (!this->m_objectiveFunction.domainSet().contains(*(this->m_initialPoint)))
215  {
216  if( m_objectiveFunction.domainSet().env().fullRank() == 0 )
217  {
218  std::cerr << "Minimization was given initial point outside of domain"
219  << std::endl;
220  }
221  queso_error();
222  }
223 
224  unsigned int dim = this->m_objectiveFunction.domainSet().vectorSpace().
225  zeroVector().sizeLocal();
226 
227  // We use m_solver_type here because we need the enum
229  {
230  this->minimize_with_gradient( dim, monitor );
231  }
232  else
233  {
234  this->minimize_no_gradient( dim, monitor );
235  }
236 
237  return;
238 }
239 
242 {
243  return this->m_objectiveFunction;
244 }
245 
246 void
248 {
249  for (unsigned int i = 0; i < initialPoint.sizeLocal(); i++) {
250  (*(this->m_initialPoint))[i] = initialPoint[i];
251  }
252 }
253 
254 const GslVector &
256 {
257  return *(this->m_minimizer);
258 }
259 
261 {
263  m_solver_type = solver;
264 }
265 
267 {
268  bool gradient_needed = false;
269 
270  switch(solver)
271  {
272  case(FLETCHER_REEVES_CG):
273  case(POLAK_RIBIERE_CG):
274  case(BFGS):
275  case(BFGS2):
276  case(STEEPEST_DESCENT):
277  {
278  gradient_needed = true;
279  break;
280  }
281  case(NELDER_MEAD):
282  case(NELDER_MEAD2):
283  case(NELDER_MEAD2_RAND):
284  {
285  break;
286  }
287  default:
288  {
289  // Wat?!
290  queso_error();
291  }
292  } // switch(solver)
293 
294  return gradient_needed;
295 }
296 
298 {
299  // Set initial point
300  gsl_vector * x = gsl_vector_alloc(dim);
301  for (unsigned int i = 0; i < dim; i++) {
302  gsl_vector_set(x, i, (*m_initialPoint)[i]);
303  }
304 
305  // Tell GSL which solver we're using
306  const gsl_multimin_fdfminimizer_type* type = NULL;
307 
308  // We use m_solver_type here because we need the enum
309  switch(m_solver_type)
310  {
311  case(FLETCHER_REEVES_CG):
312  type = gsl_multimin_fdfminimizer_conjugate_fr;
313  break;
314  case(POLAK_RIBIERE_CG):
315  type = gsl_multimin_fdfminimizer_conjugate_pr;
316  break;
317  case(BFGS):
318  type = gsl_multimin_fdfminimizer_vector_bfgs;
319  break;
320  case(BFGS2):
321  type = gsl_multimin_fdfminimizer_vector_bfgs2;
322  break;
323  case(STEEPEST_DESCENT):
324  type = gsl_multimin_fdfminimizer_steepest_descent;
325  break;
326  case(NELDER_MEAD):
327  case(NELDER_MEAD2):
328  case(NELDER_MEAD2_RAND):
329  default:
330  // Wat?!
331  queso_error();
332  }
333 
334  // Init solver
335  gsl_multimin_fdfminimizer * solver =
336  gsl_multimin_fdfminimizer_alloc(type, dim);
337 
338  // Point GSL to the right functions
339  gsl_multimin_function_fdf minusLogPosterior;
340  minusLogPosterior.n = dim;
341  minusLogPosterior.f = &c_evaluate;
342  minusLogPosterior.df = &c_evaluate_derivative;
343  minusLogPosterior.fdf = &c_evaluate_with_derivative;
344  minusLogPosterior.params = (void *)(this);
345 
346  gsl_multimin_fdfminimizer_set(solver, &minusLogPosterior, x, getFdfstepSize(), getLineTolerance());
347 
348  int status;
349  size_t iter = 0;
350 
351  do {
352  iter++;
353  status = gsl_multimin_fdfminimizer_iterate(solver);
354 
355  if (status) {
356  if( m_objectiveFunction.domainSet().env().fullRank() == 0 )
357  {
358  std::cerr << "Error while GSL does optimisation. "
359  << "See below for GSL error type." << std::endl;
360  std::cerr << "Gsl error: " << gsl_strerror(status) << std::endl;
361  }
362  break;
363  }
364 
365  status = gsl_multimin_test_gradient(solver->gradient, this->getTolerance());
366 
367  if(monitor)
368  {
369  gsl_vector* x = gsl_multimin_fdfminimizer_x(solver);
370  std::vector<double> x_min(dim);
371  for( unsigned int i = 0; i < dim; i++)
372  x_min[i] = gsl_vector_get(x,i);
373 
374  double f = gsl_multimin_fdfminimizer_minimum(solver);
375 
376  gsl_vector* grad = gsl_multimin_fdfminimizer_gradient(solver);
377  double grad_norm = gsl_blas_dnrm2(grad);
378 
379  monitor->append( x_min, f, grad_norm );
380  }
381 
382  } while ((status == GSL_CONTINUE) && (iter < this->getMaxIterations()));
383 
384  for (unsigned int i = 0; i < dim; i++) {
385  (*m_minimizer)[i] = gsl_vector_get(solver->x, i);
386  }
387 
388  // We're being good human beings and cleaning up the memory we allocated
389  gsl_multimin_fdfminimizer_free(solver);
390  gsl_vector_free(x);
391 
392  return;
393 }
394 
396 {
397  // Set initial point
398  gsl_vector* x = gsl_vector_alloc(dim);
399  for (unsigned int i = 0; i < dim; i++) {
400  gsl_vector_set(x, i, (*m_initialPoint)[i]);
401  }
402 
403  // Tell GSL which solver we're using
404  const gsl_multimin_fminimizer_type* type = NULL;
405 
406  switch(m_solver_type)
407  {
408  case(NELDER_MEAD):
409  type = gsl_multimin_fminimizer_nmsimplex;
410  break;
411  case(NELDER_MEAD2):
412  type = gsl_multimin_fminimizer_nmsimplex2;
413  break;
414  case(NELDER_MEAD2_RAND):
415  type = gsl_multimin_fminimizer_nmsimplex2rand;
416  break;
417  case(FLETCHER_REEVES_CG):
418  case(POLAK_RIBIERE_CG):
419  case(BFGS):
420  case(BFGS2):
421  case(STEEPEST_DESCENT):
422  default:
423  // Wat?!
424  queso_error();
425  }
426 
427  // Init solver
428  gsl_multimin_fminimizer* solver =
429  gsl_multimin_fminimizer_alloc(type, dim);
430 
431  // Point GSL at the right functions
432  gsl_multimin_function minusLogPosterior;
433  minusLogPosterior.n = dim;
434  minusLogPosterior.f = &c_evaluate;
435  minusLogPosterior.params = (void *)(this);
436 
437  // Needed for these gradient free algorithms.
438  gsl_vector* step_size = gsl_vector_alloc(dim);
439 
440  for(unsigned int i = 0; i < dim; i++) {
441  gsl_vector_set(step_size, i, m_fstep_size[i]);
442  }
443 
444  gsl_multimin_fminimizer_set(solver, &minusLogPosterior, x, step_size);
445 
446  int status;
447  size_t iter = 0;
448  double size = 0.0;
449 
450  do
451  {
452  iter++;
453  status = gsl_multimin_fminimizer_iterate(solver);
454 
455  if (status) {
456  if( m_objectiveFunction.domainSet().env().fullRank() == 0 )
457  {
458  std::cerr << "Error while GSL does optimisation. "
459  << "See below for GSL error type." << std::endl;
460  std::cerr << "Gsl error: " << gsl_strerror(status) << std::endl;
461  }
462  break;
463  }
464 
465  size = gsl_multimin_fminimizer_size(solver);
466 
467  status = gsl_multimin_test_size (size, this->getTolerance());
468 
469  if(monitor)
470  {
471  gsl_vector* x = gsl_multimin_fminimizer_x(solver);
472  std::vector<double> x_min(dim);
473  for( unsigned int i = 0; i < dim; i++)
474  x_min[i] = gsl_vector_get(x,i);
475 
476  double f = gsl_multimin_fminimizer_minimum(solver);
477 
478  monitor->append( x_min, f, size );
479  }
480 
481  }
482 
483  while ((status == GSL_CONTINUE) && (iter < this->getMaxIterations()));
484 
485  for (unsigned int i = 0; i < dim; i++) {
486  (*m_minimizer)[i] = gsl_vector_get(solver->x, i);
487  }
488 
489  // We're being good human beings and cleaning up the memory we allocated
490  gsl_vector_free(step_size);
491  gsl_multimin_fminimizer_free(solver);
492  gsl_vector_free(x);
493 
494  return;
495 }
496 
497 void GslOptimizer::set_step_size( const GslVector& step_size )
498 {
500  m_fstep_size = step_size;
501 }
502 
503 void GslOptimizer::set_step_size( double step_size )
504 {
506  m_fdfstep_size = step_size;
507 }
508 
510 {
511  SolverType solver_type;
512 
513  if( solver == std::string("fletcher_reeves_cg") )
514  solver_type = FLETCHER_REEVES_CG;
515  else if( solver == std::string("polak_ribiere_cg") )
516  solver_type = POLAK_RIBIERE_CG;
517  else if( solver == std::string("bfgs") )
518  solver_type = BFGS;
519  else if( solver == std::string("bfgs2") )
520  solver_type = BFGS2;
521  else if( solver == std::string("steepest_decent") ) {
523  solver_type = STEEPEST_DESCENT;
524  }
525  else if( solver == std::string("steepest_descent") )
526  solver_type = STEEPEST_DESCENT;
527  else if( solver == std::string("nelder_mead") )
528  solver_type = NELDER_MEAD;
529  else if( solver == std::string("nelder_mead2") )
530  solver_type = NELDER_MEAD2;
531  else if( solver == std::string("nelder_mead2_rand") )
532  solver_type = NELDER_MEAD2_RAND;
533  else
534  {
535  if( m_objectiveFunction.domainSet().env().fullRank() == 0 )
536  {
537  std::cerr << "Error: Invalid GslOptimizer solver name: " << solver << std::endl
538  << " Valids choices are: fletcher_reeves_cg" << std::endl
539  << " polak_ribiere_cg" << std::endl
540  << " bfgs" << std::endl
541  << " bfgs2" << std::endl
542  << " steepest_descent" << std::endl
543  << " nelder_mead" << std::endl
544  << " nelder_mead2" << std::endl
545  << " nelder_mead2_rand" << std::endl;
546  }
547  queso_error();
548  }
549 
550  return solver_type;
551 }
552 
553 void GslOptimizer::set_solver_type( std::string& solver )
554 {
556  this->set_solver_type( this->string_to_enum(solver) );
557 }
558 
559 void
560 GslOptimizer::setSolverType(std::string solverType)
561 {
562  this->m_optionsObj->m_solverType = solverType;
563  this->set_solver_type(solverType);
564 }
565 
566 void
567 GslOptimizer::setFstepSize(double fstepSize)
568 {
569  this->m_optionsObj->m_fstepSize = fstepSize;
570 
571  GslVector fstepSizeVector(
572  objectiveFunction().domainSet().vectorSpace().zeroVector());
573  fstepSizeVector.cwSet(fstepSize);
574 
575  this->set_step_size(fstepSizeVector);
576 }
577 
578 void
579 GslOptimizer::setFdfstepSize(double fdfstepSize)
580 {
581  this->m_optionsObj->m_fdfstepSize = fdfstepSize;
582  this->set_step_size(fdfstepSize);
583 }
584 
585 void
586 GslOptimizer::setLineTolerance(double lineTolerance)
587 {
588  this->m_optionsObj->m_lineTolerance = lineTolerance;
589 }
590 
591 std::string
593 {
594  return this->m_optionsObj->m_solverType;
595 }
596 
597 double
599 {
600  return this->m_optionsObj->m_fstepSize;
601 }
602 
603 double
605 {
606  return this->m_optionsObj->m_fdfstepSize;
607 }
608 
609 double
611 {
612  return this->m_optionsObj->m_lineTolerance;
613 }
614 
615 } // End namespace QUESO
virtual void setFstepSize(double fstepSize)
Sets the step size to use in gradient-free solvers.
Definition: GslOptimizer.C:567
void set_solver_type(SolverType solver)
Definition: GslOptimizer.C:260
void append(std::vector< double > &x_min, double objective, double norm)
Add current value of minimizer, objective, and error norm.
SolverType string_to_enum(std::string &solver)
Definition: GslOptimizer.C:509
A templated (base) class for handling scalar functions.
void setInitialPoint(const GslVector &intialPoint)
Set the point at which the optimization starts.
Definition: GslOptimizer.C:247
int dim
Definition: ann2fig.cpp:81
#define queso_error()
Definition: asserts.h:53
ScopedPtr< OptimizerOptions >::Type m_optionsObj
Definition: Optimizer.h:126
A base class for handling optimisation of scalar functions.
Definition: Optimizer.h:47
GslOptimizer(const BaseScalarFunction< GslVector, GslMatrix > &objectiveFunction)
Constructs an object that will maximize a scalar function.
Definition: GslOptimizer.C:159
#define queso_deprecated()
Definition: Defines.h:134
virtual void minimize(OptimizerMonitor *monitor=NULL)
Minimize the objective function, starting at m_initialPoint.
Definition: GslOptimizer.C:211
virtual void setLineTolerance(double lineTolerance)
Sets the tolerance to use for line minimisation.
Definition: GslOptimizer.C:586
Class for vector operations using GSL library.
Definition: GslVector.h:49
void minimize_with_gradient(unsigned int dim, OptimizerMonitor *monitor)
Definition: GslOptimizer.C:297
This class provides options for a Optimizer.
void c_evaluate_derivative(const gsl_vector *x, void *context, gsl_vector *derivative)
Definition: GslOptimizer.C:69
void cwSet(double value)
Component-wise sets all values to this with value.
Definition: GslVector.C:326
bool solver_needs_gradient(SolverType solver)
Helper function.
Definition: GslOptimizer.C:266
const GslVector & minimizer() const
Return the point that minimizes the objective function.
Definition: GslOptimizer.C:255
virtual ~GslOptimizer()
Destructor.
Definition: GslOptimizer.C:206
virtual double getLineTolerance() const
Gets the tolerance to use for line minimisation.
Definition: GslOptimizer.C:610
double getTolerance() const
Returns the tolerance used to test for an extremum in the optimizer.
Definition: Optimizer.C:57
const BaseScalarFunction< GslVector, GslMatrix > & m_objectiveFunction
Definition: GslOptimizer.h:155
A base class for handling optimisation of scalar functions.
Definition: GslOptimizer.h:50
SolverType m_solver_type
Definition: GslOptimizer.h:160
GslVector m_fstep_size
For use in gradient-free algorithms.
Definition: GslOptimizer.h:163
virtual void setFdfstepSize(double fdfstepSize)
Sets the step to use in gradient-based solvers.
Definition: GslOptimizer.C:579
double m_fdfstep_size
For use in gradient-based algorithms.
Definition: GslOptimizer.h:166
unsigned int sizeLocal() const
Returns the length of this vector.
Definition: GslVector.C:240
void c_evaluate_with_derivative(const gsl_vector *x, void *context, double *f, gsl_vector *derivative)
Definition: GslOptimizer.C:151
void set_step_size(const GslVector &step_size)
Sets step size used in gradient-free solvers.
Definition: GslOptimizer.C:497
void minimize_no_gradient(unsigned int dim, OptimizerMonitor *monitor)
Definition: GslOptimizer.C:395
Object to monitor convergence of optimizers.
virtual double getFdfstepSize() const
Gets the step to use in gradient-based solvers.
Definition: GslOptimizer.C:604
GslVector * m_minimizer
Definition: GslOptimizer.h:158
virtual std::string getSolverType() const
Gets the algorithm to use for minimisation.
Definition: GslOptimizer.C:592
virtual void setSolverType(std::string solverType)
Sets the algorithm to use for minimisation.
Definition: GslOptimizer.C:560
const BaseScalarFunction< GslVector, GslMatrix > & objectiveFunction() const
Returns the objective function.
Definition: GslOptimizer.C:241
double c_evaluate(const gsl_vector *x, void *context)
Definition: GslOptimizer.C:43
double getFiniteDifferenceStepSize() const
Returns the step size used in the finite difference formula.
Definition: Optimizer.C:63
virtual double getFstepSize() const
Gets the step size to use in gradient-free solvers.
Definition: GslOptimizer.C:598
unsigned int getMaxIterations() const
Returns the maximum number of iterations the optimizer will do.
Definition: Optimizer.C:51
GslVector * m_initialPoint
Definition: GslOptimizer.h:157

Generated on Thu Dec 15 2016 13:23:10 for queso-0.56.1 by  doxygen 1.8.5