queso-0.53.0
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 <gsl/gsl_multimin.h>
28 #include <gsl/gsl_blas.h>
29 #include <queso/GslVector.h>
30 #include <queso/VectorSpace.h>
31 #include <queso/ScalarFunction.h>
32 #include <queso/GslOptimizer.h>
33 #include <queso/OptimizerMonitor.h>
34 
35 namespace QUESO {
36 
37 // We need to extern "C" because gsl needs a pointer to a C function to
38 // minimize
39 extern "C" {
40  // This evaluate -log posterior
41  double c_evaluate(const gsl_vector * x, void * context) {
42 
43  GslOptimizer * optimizer = static_cast<GslOptimizer * >(context);
44 
45  GslVector state(
46  optimizer->objectiveFunction().domainSet().vectorSpace().zeroVector());
47 
48  // DM: Doing this copy sucks, but whatever. It'll do for now.
49  for (unsigned int i = 0; i < state.sizeLocal(); i++) {
50  state[i] = gsl_vector_get(x, i);
51  }
52 
53  // Bail early if GSL tries to evaluate outside of the domain
54  if (!optimizer->objectiveFunction().domainSet().contains(state)) {
55  return GSL_NAN;
56  }
57 
58  // Should cache derivative here so we don't a) segfault in the user's code
59  // and b) so we don't recompute stuff
60  double result = -optimizer->objectiveFunction().lnValue(state, NULL, NULL,
61  NULL, NULL);
62 
63  return result;
64  }
65 
66  // This evaluates the derivative of -log posterior
67  void c_evaluate_derivative(const gsl_vector * x, void * context,
68  gsl_vector * derivative) {
69  GslOptimizer * optimizer = static_cast<GslOptimizer * >(context);
70 
71  GslVector state(
72  optimizer->objectiveFunction().domainSet().vectorSpace().zeroVector());
73  GslVector deriv(
74  optimizer->objectiveFunction().domainSet().vectorSpace().zeroVector());
75 
76  // DM: Doing this copy sucks, but whatever. It'll do for now.
77  for (unsigned int i = 0; i < state.sizeLocal(); i++) {
78  state[i] = gsl_vector_get(x, i);
79 
80  // We fill with GSL_NAN and use it as a flag to check later that the user
81  // actually fills the derivative vector with stuff
82  deriv[i] = GSL_NAN;
83  }
84 
85  if (!optimizer->objectiveFunction().domainSet().contains(state)) {
86  // Fill derivative with error codes if the point is outside of the
87  // domain
88  for (unsigned int i = 0; i < deriv.sizeLocal(); i++) {
89  gsl_vector_set(derivative, i, GSL_NAN);
90  }
91  }
92  else {
93  // We should cache the return value here so we don't recompute stuff
94  double fx = -optimizer->objectiveFunction().lnValue(state, NULL, &deriv,
95  NULL, NULL);
96 
97  // Decide whether or not we need to do a finite difference based on
98  // whether the user actually filled deriv with values that are not
99  // GSL_NAN
100  //
101  // We're currently doing this check every time this function gets called.
102  // We could probably pull this logic out of here and put it somewhere
103  // where it only happens once
104  bool userComputedDerivative = true;
105  for (unsigned int i = 0; i < deriv.sizeLocal(); i++) {
106  // If the user missed out a derivative in any direction, fall back to
107  // a finite difference
108  if (gsl_isnan(deriv[i])) {
109  userComputedDerivative = false;
110  break;
111  }
112  }
113 
114  if (userComputedDerivative) {
115  for (unsigned int i = 0; i < deriv.sizeLocal(); i++) {
116  gsl_vector_set(derivative, i, -deriv[i]); // We need the minus sign
117  }
118  }
119  else {
120  // Finite difference step-size
121  double h = optimizer->getFiniteDifferenceStepSize();
122 
123  // User did not provide a derivative, so do a finite difference
124  for (unsigned int i = 0; i < deriv.sizeLocal(); i++) {
125  double tempState = state[i];
126  state[i] += h;
127 
128  // User didn't provide a derivative, so we don't bother passing in
129  // the derivative vector again
130  double fxph = -optimizer->objectiveFunction().lnValue(state, NULL,
131  NULL, NULL, NULL);
132 
133  // Reset the state back to what it was before
134  state[i] = tempState;
135 
136  // Make sure we didn't do anything dumb and tell gsl if we did
137  if (!gsl_isnan(fx) && !gsl_isnan(fxph)) {
138  gsl_vector_set(derivative, i, (fxph - fx) / h);
139  }
140  else {
141  gsl_vector_set(derivative, i, GSL_NAN);
142  }
143  }
144  }
145  }
146  }
147 
148  // This evaluates -log posterior and the derivative of -log posterior
149  void c_evaluate_with_derivative(const gsl_vector * x, void * context,
150  double * f, gsl_vector * derivative) {
151  // We don't need to call both of these
152  *f = c_evaluate(x, context);
153  c_evaluate_derivative(x, context, derivative);
154  }
155 } // End extern "C"
156 
158  const BaseScalarFunction<GslVector, GslMatrix> & objectiveFunction)
159  : BaseOptimizer(),
160  m_objectiveFunction(objectiveFunction),
161  m_initialPoint(new GslVector(objectiveFunction.domainSet().
162  vectorSpace().zeroVector())),
163  m_minimizer(new GslVector(this->m_objectiveFunction.domainSet().
164  vectorSpace().zeroVector())),
165  m_solver_type(BFGS2),
166  m_fstep_size(this->m_objectiveFunction.domainSet().vectorSpace().zeroVector()),
167  m_fdfstep_size(1.0),
168  m_line_tol(0.1)
169 {
170  // We initialize the minimizer to GSL_NAN just in case the optimization fails
171  m_minimizer->cwSet(GSL_NAN);
172 
173  // Set to documented default value.
174  m_fstep_size.cwSet(0.1);
175 }
176 
178 {
179  delete this->m_initialPoint;
180 }
181 
183 
184  // First check that initial guess is reasonable
185  if (!this->m_objectiveFunction.domainSet().contains(*(this->m_initialPoint)))
186  {
187  if( m_objectiveFunction.domainSet().env().fullRank() == 0 )
188  {
189  std::cerr << "Minimization was given initial point outside of domain"
190  << std::endl;
191  }
192  queso_error();
193  }
194 
195  unsigned int dim = this->m_objectiveFunction.domainSet().vectorSpace().
196  zeroVector().sizeLocal();
197 
199  {
200  this->minimize_with_gradient( dim, monitor );
201  }
202  else
203  {
204  this->minimize_no_gradient( dim, monitor );
205  }
206 
207  return;
208 }
209 
212 {
213  return this->m_objectiveFunction;
214 }
215 
216 void
218 {
219  for (unsigned int i = 0; i < initialPoint.sizeLocal(); i++) {
220  (*(this->m_initialPoint))[i] = initialPoint[i];
221  }
222 }
223 
224 const GslVector &
226 {
227  return *(this->m_minimizer);
228 }
229 
231  {
232  m_solver_type = solver;
233  }
234 
236  {
237  bool gradient_needed = false;
238 
239  switch(solver)
240  {
241  case(FLETCHER_REEVES_CG):
242  case(POLAK_RIBIERE_CG):
243  case(BFGS):
244  case(BFGS2):
245  case(STEEPEST_DESCENT):
246  {
247  gradient_needed = true;
248  break;
249  }
250  case(NELDER_MEAD):
251  case(NELDER_MEAD2):
252  case(NELDER_MEAD2_RAND):
253  {
254  break;
255  }
256  default:
257  {
258  // Wat?!
259  queso_error();
260  }
261  } // switch(solver)
262 
263  return gradient_needed;
264  }
265 
267  {
268  // Set initial point
269  gsl_vector * x = gsl_vector_alloc(dim);
270  for (unsigned int i = 0; i < dim; i++) {
271  gsl_vector_set(x, i, (*m_initialPoint)[i]);
272  }
273 
274  // Tell GSL which solver we're using
275  const gsl_multimin_fdfminimizer_type* type = NULL;
276 
277  switch(m_solver_type)
278  {
279  case(FLETCHER_REEVES_CG):
280  type = gsl_multimin_fdfminimizer_conjugate_fr;
281  break;
282  case(POLAK_RIBIERE_CG):
283  type = gsl_multimin_fdfminimizer_conjugate_pr;
284  break;
285  case(BFGS):
286  type = gsl_multimin_fdfminimizer_vector_bfgs;
287  break;
288  case(BFGS2):
289  type = gsl_multimin_fdfminimizer_vector_bfgs2;
290  break;
291  case(STEEPEST_DESCENT):
292  type = gsl_multimin_fdfminimizer_steepest_descent;
293  break;
294  case(NELDER_MEAD):
295  case(NELDER_MEAD2):
296  case(NELDER_MEAD2_RAND):
297  default:
298  // Wat?!
299  queso_error();
300  }
301 
302  // Init solver
303  gsl_multimin_fdfminimizer * solver =
304  gsl_multimin_fdfminimizer_alloc(type, dim);
305 
306  // Point GSL to the right functions
307  gsl_multimin_function_fdf minusLogPosterior;
308  minusLogPosterior.n = dim;
309  minusLogPosterior.f = &c_evaluate;
310  minusLogPosterior.df = &c_evaluate_derivative;
311  minusLogPosterior.fdf = &c_evaluate_with_derivative;
312  minusLogPosterior.params = (void *)(this);
313 
314  gsl_multimin_fdfminimizer_set(solver, &minusLogPosterior, x, m_fdfstep_size, m_line_tol);
315 
316  int status;
317  size_t iter = 0;
318 
319  do {
320  iter++;
321  status = gsl_multimin_fdfminimizer_iterate(solver);
322 
323  if (status) {
324  if( m_objectiveFunction.domainSet().env().fullRank() == 0 )
325  {
326  std::cerr << "Error while GSL does optimisation. "
327  << "See below for GSL error type." << std::endl;
328  std::cerr << "Gsl error: " << gsl_strerror(status) << std::endl;
329  }
330  break;
331  }
332 
333  status = gsl_multimin_test_gradient(solver->gradient, this->getTolerance());
334 
335  if(monitor)
336  {
337  gsl_vector* x = gsl_multimin_fdfminimizer_x(solver);
338  std::vector<double> x_min(dim);
339  for( unsigned int i = 0; i < dim; i++)
340  x_min[i] = gsl_vector_get(x,i);
341 
342  double f = gsl_multimin_fdfminimizer_minimum(solver);
343 
344  gsl_vector* grad = gsl_multimin_fdfminimizer_gradient(solver);
345  double grad_norm = gsl_blas_dnrm2(grad);
346 
347  monitor->append( x_min, f, grad_norm );
348  }
349 
350  } while ((status == GSL_CONTINUE) && (iter < this->getMaxIterations()));
351 
352  for (unsigned int i = 0; i < dim; i++) {
353  (*m_minimizer)[i] = gsl_vector_get(solver->x, i);
354  }
355 
356  // We're being good human beings and cleaning up the memory we allocated
357  gsl_multimin_fdfminimizer_free(solver);
358  gsl_vector_free(x);
359 
360  return;
361  }
362 
364  {
365  // Set initial point
366  gsl_vector* x = gsl_vector_alloc(dim);
367  for (unsigned int i = 0; i < dim; i++) {
368  gsl_vector_set(x, i, (*m_initialPoint)[i]);
369  }
370 
371  // Tell GSL which solver we're using
372  const gsl_multimin_fminimizer_type* type = NULL;
373 
374  switch(m_solver_type)
375  {
376  case(NELDER_MEAD):
377  type = gsl_multimin_fminimizer_nmsimplex;
378  break;
379  case(NELDER_MEAD2):
380  type = gsl_multimin_fminimizer_nmsimplex2;
381  break;
382  case(NELDER_MEAD2_RAND):
383  type = gsl_multimin_fminimizer_nmsimplex2rand;
384  break;
385  case(FLETCHER_REEVES_CG):
386  case(POLAK_RIBIERE_CG):
387  case(BFGS):
388  case(BFGS2):
389  case(STEEPEST_DESCENT):
390  default:
391  // Wat?!
392  queso_error();
393  }
394 
395  // Init solver
396  gsl_multimin_fminimizer* solver =
397  gsl_multimin_fminimizer_alloc(type, dim);
398 
399  // Point GSL at the right functions
400  gsl_multimin_function minusLogPosterior;
401  minusLogPosterior.n = dim;
402  minusLogPosterior.f = &c_evaluate;
403  minusLogPosterior.params = (void *)(this);
404 
405  // Needed for these gradient free algorithms.
406  gsl_vector* step_size = gsl_vector_alloc(dim);
407 
408  for(unsigned int i = 0; i < dim; i++) {
409  gsl_vector_set(step_size, i, m_fstep_size[i]);
410  }
411 
412  gsl_multimin_fminimizer_set(solver, &minusLogPosterior, x, step_size);
413 
414  int status;
415  size_t iter = 0;
416  double size = 0.0;
417 
418  do
419  {
420  iter++;
421  status = gsl_multimin_fminimizer_iterate(solver);
422 
423  if (status) {
424  if( m_objectiveFunction.domainSet().env().fullRank() == 0 )
425  {
426  std::cerr << "Error while GSL does optimisation. "
427  << "See below for GSL error type." << std::endl;
428  std::cerr << "Gsl error: " << gsl_strerror(status) << std::endl;
429  }
430  break;
431  }
432 
433  size = gsl_multimin_fminimizer_size(solver);
434 
435  status = gsl_multimin_test_size (size, this->getTolerance());
436 
437  if(monitor)
438  {
439  gsl_vector* x = gsl_multimin_fminimizer_x(solver);
440  std::vector<double> x_min(dim);
441  for( unsigned int i = 0; i < dim; i++)
442  x_min[i] = gsl_vector_get(x,i);
443 
444  double f = gsl_multimin_fminimizer_minimum(solver);
445 
446  monitor->append( x_min, f, size );
447  }
448 
449  }
450 
451  while ((status == GSL_CONTINUE) && (iter < this->getMaxIterations()));
452 
453  for (unsigned int i = 0; i < dim; i++) {
454  (*m_minimizer)[i] = gsl_vector_get(solver->x, i);
455  }
456 
457  // We're being good human beings and cleaning up the memory we allocated
458  gsl_vector_free(step_size);
459  gsl_multimin_fminimizer_free(solver);
460  gsl_vector_free(x);
461 
462  return;
463  }
464 
465  void GslOptimizer::set_step_size( const GslVector& step_size )
466  {
467  m_fstep_size = step_size;
468  }
469 
470  void GslOptimizer::set_step_size( double step_size )
471  {
472  m_fdfstep_size = step_size;
473  }
474 
476  {
477  SolverType solver_type;
478 
479  if( solver == std::string("fletcher_reeves_cg") )
480  solver_type = FLETCHER_REEVES_CG;
481  else if( solver == std::string("polak_ribiere_cg") )
482  solver_type = POLAK_RIBIERE_CG;
483  else if( solver == std::string("bfgs") )
484  solver_type = BFGS;
485  else if( solver == std::string("bfgs2") )
486  solver_type = BFGS2;
487  else if( solver == std::string("steepest_decent") )
488  solver_type = STEEPEST_DESCENT;
489  else if( solver == std::string("nelder_mead") )
490  solver_type = NELDER_MEAD;
491  else if( solver == std::string("nelder_mead2") )
492  solver_type = NELDER_MEAD2;
493  else if( solver == std::string("nelder_mead2_rand") )
494  solver_type = NELDER_MEAD2_RAND;
495  else
496  {
497  if( m_objectiveFunction.domainSet().env().fullRank() == 0 )
498  {
499  std::cerr << "Error: Invalid GslOptimizer solver name: " << solver << std::endl
500  << " Valids choices are: fletcher_reeves_cg" << std::endl
501  << " polak_ribiere_cg" << std::endl
502  << " bfgs" << std::endl
503  << " bfgs2" << std::endl
504  << " steepest_decent" << std::endl
505  << " nelder_mead" << std::endl
506  << " nelder_mead2" << std::endl
507  << " nelder_mead2_rand" << std::endl;
508  }
509  queso_error();
510  }
511 
512  return solver_type;
513  }
514 
515  void GslOptimizer::set_solver_type( std::string& solver )
516  {
517  this->set_solver_type( this->string_to_enum(solver) );
518  }
519 
520 } // End namespace QUESO
GslVector * m_minimizer
Definition: GslOptimizer.h:123
void set_solver_type(SolverType solver)
Definition: GslOptimizer.C:230
void setInitialPoint(const GslVector &intialPoint)
Set the point at which the optimization starts.
Definition: GslOptimizer.C:217
virtual ~GslOptimizer()
Destructor.
Definition: GslOptimizer.C:177
GslOptimizer(const BaseScalarFunction< GslVector, GslMatrix > &objectiveFunction)
Constructs an object that will maximize a scalar function.
Definition: GslOptimizer.C:157
void set_step_size(const GslVector &step_size)
Sets step size used in gradient-free solvers.
Definition: GslOptimizer.C:465
double m_line_tol
Line minimization tolerance in gradient-based algorithms.
Definition: GslOptimizer.h:134
GslVector m_fstep_size
For use in gradient-free algorithms.
Definition: GslOptimizer.h:128
void append(std::vector< double > &x_min, double objective, double norm)
Add current value of minimizer, objective, and error norm.
unsigned int sizeLocal() const
Returns the length of this vector.
Definition: GslVector.C:253
A templated (base) class for handling scalar functions.
SolverType string_to_enum(std::string &solver)
Definition: GslOptimizer.C:475
A base class for handling optimisation of scalar functions.
Definition: Optimizer.h:44
void c_evaluate_derivative(const gsl_vector *x, void *context, gsl_vector *derivative)
Definition: GslOptimizer.C:67
GslVector * m_initialPoint
Definition: GslOptimizer.h:122
virtual void minimize(OptimizerMonitor *monitor=NULL)
Minimize the objective function, starting at m_initialPoint.
Definition: GslOptimizer.C:182
double c_evaluate(const gsl_vector *x, void *context)
Definition: GslOptimizer.C:41
#define queso_error()
Definition: asserts.h:53
unsigned int getMaxIterations() const
Returns the maximum number of iterations the optimizer will do.
Definition: Optimizer.C:41
Class for vector operations using GSL library.
Definition: GslVector.h:48
void c_evaluate_with_derivative(const gsl_vector *x, void *context, double *f, gsl_vector *derivative)
Definition: GslOptimizer.C:149
void cwSet(double value)
Component-wise sets all values to this with value.
Definition: GslVector.C:339
double getTolerance() const
Returns the tolerance used to test for an extremum in the optimizer.
Definition: Optimizer.C:47
bool solver_needs_gradient(SolverType solver)
Helper function.
Definition: GslOptimizer.C:235
Object to monitor convergence of optimizers.
const BaseScalarFunction< GslVector, GslMatrix > & m_objectiveFunction
Definition: GslOptimizer.h:120
const BaseScalarFunction< GslVector, GslMatrix > & objectiveFunction() const
Returns the objective function.
Definition: GslOptimizer.C:211
void minimize_no_gradient(unsigned int dim, OptimizerMonitor *monitor)
Definition: GslOptimizer.C:363
SolverType m_solver_type
Definition: GslOptimizer.h:125
double m_fdfstep_size
For use in gradient-based algorithms.
Definition: GslOptimizer.h:131
double getFiniteDifferenceStepSize() const
Returns the step size used in the finite difference formula.
Definition: Optimizer.C:53
void minimize_with_gradient(unsigned int dim, OptimizerMonitor *monitor)
Definition: GslOptimizer.C:266
int dim
Definition: ann2fig.cpp:81
A base class for handling optimisation of scalar functions.
Definition: GslOptimizer.h:50
const GslVector & minimizer() const
Return the point that minimizes the objective function.
Definition: GslOptimizer.C:225

Generated on Thu Jun 11 2015 13:52:31 for queso-0.53.0 by  doxygen 1.8.5