queso-0.57.0
Public Types | Public Member Functions | Private Member Functions | Private Attributes | List of all members
QUESO::GslOptimizer Class Reference

A base class for handling optimisation of scalar functions. More...

#include <GslOptimizer.h>

Inheritance diagram for QUESO::GslOptimizer:
Inheritance graph
[legend]
Collaboration diagram for QUESO::GslOptimizer:
Collaboration graph
[legend]

Public Types

enum  SolverType {
  FLETCHER_REEVES_CG, POLAK_RIBIERE_CG, BFGS, BFGS2,
  STEEPEST_DESCENT, NELDER_MEAD, NELDER_MEAD2, NELDER_MEAD2_RAND
}
 

Public Member Functions

 GslOptimizer (const BaseScalarFunction< GslVector, GslMatrix > &objectiveFunction)
 Constructs an object that will maximize a scalar function. More...
 
 GslOptimizer (OptimizerOptions options, const BaseScalarFunction< GslVector, GslMatrix > &objectiveFunction)
 Constructs an object that will maximize a scalar function. More...
 
virtual ~GslOptimizer ()
 Destructor. More...
 
virtual void minimize (OptimizerMonitor *monitor=NULL)
 Minimize the objective function, starting at m_initialPoint. More...
 
const BaseScalarFunction
< GslVector, GslMatrix > & 
objectiveFunction () const
 Returns the objective function. More...
 
void setInitialPoint (const GslVector &intialPoint)
 Set the point at which the optimization starts. More...
 
const GslVectorminimizer () const
 Return the point that minimizes the objective function. More...
 
void set_solver_type (SolverType solver)
 
void set_solver_type (std::string &solver)
 
SolverType string_to_enum (std::string &solver)
 
void set_step_size (const GslVector &step_size)
 Sets step size used in gradient-free solvers. More...
 
void set_step_size (double step_size)
 Sets step size used in gradient-based solvers. More...
 
void set_line_tol (double tol)
 Set GSL line minimization tolerance. More...
 
virtual std::string getSolverType () const
 Gets the algorithm to use for minimisation. More...
 
virtual double getFstepSize () const
 Gets the step size to use in gradient-free solvers. More...
 
virtual double getFdfstepSize () const
 Gets the step to use in gradient-based solvers. More...
 
virtual double getLineTolerance () const
 Gets the tolerance to use for line minimisation. More...
 
virtual void setSolverType (std::string solverType)
 Sets the algorithm to use for minimisation. More...
 
virtual void setFstepSize (double fstepSize)
 Sets the step size to use in gradient-free solvers. More...
 
virtual void setFdfstepSize (double fdfstepSize)
 Sets the step to use in gradient-based solvers. More...
 
virtual void setLineTolerance (double lineTolerance)
 Sets the tolerance to use for line minimisation. More...
 
- Public Member Functions inherited from QUESO::BaseOptimizer
 BaseOptimizer ()
 Default constructor. More...
 
 BaseOptimizer (OptimizerOptions options)
 Constructor that takes an options object. More...
 
virtual ~BaseOptimizer ()
 Destructor. More...
 
unsigned int getMaxIterations () const
 Returns the maximum number of iterations the optimizer will do. More...
 
double getTolerance () const
 Returns the tolerance used to test for an extremum in the optimizer. More...
 
double getFiniteDifferenceStepSize () const
 Returns the step size used in the finite difference formula. More...
 
void setMaxIterations (unsigned int maxIterations)
 Sets the maximum number of iterations to be used by the optimizer. More...
 
void setTolerance (double tolerance)
 Sets the tolerance the optimizer will use to test for an extremum. More...
 
void setFiniteDifferenceStepSize (double h)
 Sets the step to use in the finite difference derivative. More...
 

Private Member Functions

bool solver_needs_gradient (SolverType solver)
 Helper function. More...
 
void minimize_with_gradient (unsigned int dim, OptimizerMonitor *monitor)
 
void minimize_no_gradient (unsigned int dim, OptimizerMonitor *monitor)
 

Private Attributes

const BaseScalarFunction
< GslVector, GslMatrix > & 
m_objectiveFunction
 
GslVectorm_initialPoint
 
GslVectorm_minimizer
 
SolverType m_solver_type
 
GslVector m_fstep_size
 For use in gradient-free algorithms. More...
 
double m_fdfstep_size
 For use in gradient-based algorithms. More...
 
double m_line_tol
 Line minimization tolerance in gradient-based algorithms. More...
 

Additional Inherited Members

- Protected Attributes inherited from QUESO::BaseOptimizer
unsigned int m_maxIterations
 
double m_tolerance
 
double m_finiteDifferenceStepSize
 
std::string m_solverType
 
double m_fstepSize
 
double m_fdfstepSize
 
double m_lineTolerance
 
ScopedPtr< OptimizerOptions >::Type m_optionsObj
 

Detailed Description

A base class for handling optimisation of scalar functions.

WRITE DOCS HERE

Definition at line 50 of file GslOptimizer.h.

Member Enumeration Documentation

Constructor & Destructor Documentation

QUESO::GslOptimizer::GslOptimizer ( const BaseScalarFunction< GslVector, GslMatrix > &  objectiveFunction)

Constructs an object that will maximize a scalar function.

The function objectiveFunction is the function that will be maximized.

Definition at line 155 of file GslOptimizer.C.

References QUESO::GslVector::cwSet(), getFstepSize(), getSolverType(), m_fstep_size, m_minimizer, and setSolverType().

157  : BaseOptimizer(),
159  m_initialPoint(new GslVector(objectiveFunction.domainSet().
160  vectorSpace().zeroVector())),
161  m_minimizer(new GslVector(this->m_objectiveFunction.domainSet().
162  vectorSpace().zeroVector())),
164  m_fstep_size(this->m_objectiveFunction.domainSet().vectorSpace().zeroVector()),
167 {
168  // We initialize the minimizer to GSL_NAN just in case the optimization fails
169  m_minimizer->cwSet(GSL_NAN);
170 
171  // Set to documented default value.
173 
174  // Set solver type to the one set in the options object
176 }
BaseOptimizer()
Default constructor.
Definition: Optimizer.C:29
GslVector * m_initialPoint
Definition: GslOptimizer.h:157
virtual double getLineTolerance() const
Gets the tolerance to use for line minimisation.
Definition: GslOptimizer.C:606
const BaseScalarFunction< GslVector, GslMatrix > & objectiveFunction() const
Returns the objective function.
Definition: GslOptimizer.C:237
void cwSet(double value)
Component-wise sets all values to this with value.
Definition: GslVector.C:327
SolverType m_solver_type
Definition: GslOptimizer.h:160
double m_fdfstep_size
For use in gradient-based algorithms.
Definition: GslOptimizer.h:166
const BaseScalarFunction< GslVector, GslMatrix > & m_objectiveFunction
Definition: GslOptimizer.h:155
GslVector * m_minimizer
Definition: GslOptimizer.h:158
double m_line_tol
Line minimization tolerance in gradient-based algorithms.
Definition: GslOptimizer.h:169
virtual std::string getSolverType() const
Gets the algorithm to use for minimisation.
Definition: GslOptimizer.C:588
virtual void setSolverType(std::string solverType)
Sets the algorithm to use for minimisation.
Definition: GslOptimizer.C:556
GslVector m_fstep_size
For use in gradient-free algorithms.
Definition: GslOptimizer.h:163
virtual double getFstepSize() const
Gets the step size to use in gradient-free solvers.
Definition: GslOptimizer.C:594
virtual double getFdfstepSize() const
Gets the step to use in gradient-based solvers.
Definition: GslOptimizer.C:600
QUESO::GslOptimizer::GslOptimizer ( OptimizerOptions  options,
const BaseScalarFunction< GslVector, GslMatrix > &  objectiveFunction 
)

Constructs an object that will maximize a scalar function.

The function objectiveFunction is the function that will be maximized. This constructor allows the passing of custom options to optimizer to modify things like tolerance, maximum number of iterations, and finite difference step size.

Definition at line 178 of file GslOptimizer.C.

References QUESO::GslVector::cwSet(), getFstepSize(), getSolverType(), m_fstep_size, m_minimizer, and setSolverType().

181  : BaseOptimizer(options),
183  m_initialPoint(new GslVector(objectiveFunction.domainSet().
184  vectorSpace().zeroVector())),
185  m_minimizer(new GslVector(this->m_objectiveFunction.domainSet().
186  vectorSpace().zeroVector())),
188  m_fstep_size(this->m_objectiveFunction.domainSet().vectorSpace().zeroVector()),
191 {
192  // We initialize the minimizer to GSL_NAN just in case the optimization fails
193  m_minimizer->cwSet(GSL_NAN);
194 
195  // Set to documented default value.
197 
198  // Set solver type to the one set in the options object
200 }
BaseOptimizer()
Default constructor.
Definition: Optimizer.C:29
GslVector * m_initialPoint
Definition: GslOptimizer.h:157
virtual double getLineTolerance() const
Gets the tolerance to use for line minimisation.
Definition: GslOptimizer.C:606
const BaseScalarFunction< GslVector, GslMatrix > & objectiveFunction() const
Returns the objective function.
Definition: GslOptimizer.C:237
void cwSet(double value)
Component-wise sets all values to this with value.
Definition: GslVector.C:327
SolverType m_solver_type
Definition: GslOptimizer.h:160
double m_fdfstep_size
For use in gradient-based algorithms.
Definition: GslOptimizer.h:166
const BaseScalarFunction< GslVector, GslMatrix > & m_objectiveFunction
Definition: GslOptimizer.h:155
GslVector * m_minimizer
Definition: GslOptimizer.h:158
double m_line_tol
Line minimization tolerance in gradient-based algorithms.
Definition: GslOptimizer.h:169
virtual std::string getSolverType() const
Gets the algorithm to use for minimisation.
Definition: GslOptimizer.C:588
virtual void setSolverType(std::string solverType)
Sets the algorithm to use for minimisation.
Definition: GslOptimizer.C:556
GslVector m_fstep_size
For use in gradient-free algorithms.
Definition: GslOptimizer.h:163
virtual double getFstepSize() const
Gets the step size to use in gradient-free solvers.
Definition: GslOptimizer.C:594
virtual double getFdfstepSize() const
Gets the step to use in gradient-based solvers.
Definition: GslOptimizer.C:600
QUESO::GslOptimizer::~GslOptimizer ( )
virtual

Destructor.

Definition at line 202 of file GslOptimizer.C.

References m_initialPoint.

203 {
204  delete this->m_initialPoint;
205 }
GslVector * m_initialPoint
Definition: GslOptimizer.h:157

Member Function Documentation

double QUESO::GslOptimizer::getFdfstepSize ( ) const
virtual

Gets the step to use in gradient-based solvers.

Reimplemented from QUESO::BaseOptimizer.

Definition at line 600 of file GslOptimizer.C.

References QUESO::BaseOptimizer::m_optionsObj.

Referenced by minimize_with_gradient().

601 {
602  return this->m_optionsObj->m_fdfstepSize;
603 }
ScopedPtr< OptimizerOptions >::Type m_optionsObj
Definition: Optimizer.h:126
double QUESO::GslOptimizer::getFstepSize ( ) const
virtual

Gets the step size to use in gradient-free solvers.

Reimplemented from QUESO::BaseOptimizer.

Definition at line 594 of file GslOptimizer.C.

References QUESO::BaseOptimizer::m_optionsObj.

Referenced by GslOptimizer().

595 {
596  return this->m_optionsObj->m_fstepSize;
597 }
ScopedPtr< OptimizerOptions >::Type m_optionsObj
Definition: Optimizer.h:126
double QUESO::GslOptimizer::getLineTolerance ( ) const
virtual

Gets the tolerance to use for line minimisation.

Reimplemented from QUESO::BaseOptimizer.

Definition at line 606 of file GslOptimizer.C.

References QUESO::BaseOptimizer::m_optionsObj.

Referenced by minimize_with_gradient().

607 {
608  return this->m_optionsObj->m_lineTolerance;
609 }
ScopedPtr< OptimizerOptions >::Type m_optionsObj
Definition: Optimizer.h:126
std::string QUESO::GslOptimizer::getSolverType ( ) const
virtual

Gets the algorithm to use for minimisation.

Reimplemented from QUESO::BaseOptimizer.

Definition at line 588 of file GslOptimizer.C.

References QUESO::BaseOptimizer::m_optionsObj.

Referenced by GslOptimizer().

589 {
590  return this->m_optionsObj->m_solverType;
591 }
ScopedPtr< OptimizerOptions >::Type m_optionsObj
Definition: Optimizer.h:126
void QUESO::GslOptimizer::minimize ( OptimizerMonitor monitor = NULL)
virtual

Minimize the objective function, starting at m_initialPoint.

m_initialPoint is handled in the derived class

Implements QUESO::BaseOptimizer.

Definition at line 207 of file GslOptimizer.C.

References dim, m_initialPoint, m_objectiveFunction, m_solver_type, minimize_no_gradient(), minimize_with_gradient(), and solver_needs_gradient().

Referenced by QUESO::StatisticalInverseProblem< P_V, P_M >::solveWithBayesMetropolisHastings().

207  {
208 
209  // First check that initial guess is reasonable
210  if (!this->m_objectiveFunction.domainSet().contains(*(this->m_initialPoint)))
211  {
212  if( m_objectiveFunction.domainSet().env().fullRank() == 0 )
213  {
214  std::cerr << "Minimization was given initial point outside of domain"
215  << std::endl;
216  }
217  queso_error();
218  }
219 
220  unsigned int dim = this->m_objectiveFunction.domainSet().vectorSpace().
221  zeroVector().sizeLocal();
222 
223  // We use m_solver_type here because we need the enum
225  {
226  this->minimize_with_gradient( dim, monitor );
227  }
228  else
229  {
230  this->minimize_no_gradient( dim, monitor );
231  }
232 
233  return;
234 }
bool solver_needs_gradient(SolverType solver)
Helper function.
Definition: GslOptimizer.C:262
GslVector * m_initialPoint
Definition: GslOptimizer.h:157
SolverType m_solver_type
Definition: GslOptimizer.h:160
const BaseScalarFunction< GslVector, GslMatrix > & m_objectiveFunction
Definition: GslOptimizer.h:155
void minimize_no_gradient(unsigned int dim, OptimizerMonitor *monitor)
Definition: GslOptimizer.C:391
void minimize_with_gradient(unsigned int dim, OptimizerMonitor *monitor)
Definition: GslOptimizer.C:293
int dim
Definition: ann2fig.cpp:81
void QUESO::GslOptimizer::minimize_no_gradient ( unsigned int  dim,
OptimizerMonitor monitor 
)
private

Definition at line 391 of file GslOptimizer.C.

References QUESO::OptimizerMonitor::append(), BFGS, BFGS2, QUESO::c_evaluate(), dim, FLETCHER_REEVES_CG, QUESO::BaseOptimizer::getMaxIterations(), QUESO::BaseOptimizer::getTolerance(), m_fstep_size, m_initialPoint, m_objectiveFunction, m_solver_type, NELDER_MEAD, NELDER_MEAD2, NELDER_MEAD2_RAND, POLAK_RIBIERE_CG, QUESO::size, and STEEPEST_DESCENT.

Referenced by minimize().

392 {
393  // Set initial point
394  gsl_vector* x = gsl_vector_alloc(dim);
395  for (unsigned int i = 0; i < dim; i++) {
396  gsl_vector_set(x, i, (*m_initialPoint)[i]);
397  }
398 
399  // Tell GSL which solver we're using
400  const gsl_multimin_fminimizer_type* type = NULL;
401 
402  switch(m_solver_type)
403  {
404  case(NELDER_MEAD):
405  type = gsl_multimin_fminimizer_nmsimplex;
406  break;
407  case(NELDER_MEAD2):
408  type = gsl_multimin_fminimizer_nmsimplex2;
409  break;
410  case(NELDER_MEAD2_RAND):
411  type = gsl_multimin_fminimizer_nmsimplex2rand;
412  break;
413  case(FLETCHER_REEVES_CG):
414  case(POLAK_RIBIERE_CG):
415  case(BFGS):
416  case(BFGS2):
417  case(STEEPEST_DESCENT):
418  default:
419  // Wat?!
420  queso_error();
421  }
422 
423  // Init solver
424  gsl_multimin_fminimizer* solver =
425  gsl_multimin_fminimizer_alloc(type, dim);
426 
427  // Point GSL at the right functions
428  gsl_multimin_function minusLogPosterior;
429  minusLogPosterior.n = dim;
430  minusLogPosterior.f = &c_evaluate;
431  minusLogPosterior.params = (void *)(this);
432 
433  // Needed for these gradient free algorithms.
434  gsl_vector* step_size = gsl_vector_alloc(dim);
435 
436  for(unsigned int i = 0; i < dim; i++) {
437  gsl_vector_set(step_size, i, m_fstep_size[i]);
438  }
439 
440  gsl_multimin_fminimizer_set(solver, &minusLogPosterior, x, step_size);
441 
442  int status;
443  size_t iter = 0;
444  double size = 0.0;
445 
446  do
447  {
448  iter++;
449  status = gsl_multimin_fminimizer_iterate(solver);
450 
451  if (status) {
452  if( m_objectiveFunction.domainSet().env().fullRank() == 0 )
453  {
454  std::cerr << "Error while GSL does optimisation. "
455  << "See below for GSL error type." << std::endl;
456  std::cerr << "Gsl error: " << gsl_strerror(status) << std::endl;
457  }
458  break;
459  }
460 
461  size = gsl_multimin_fminimizer_size(solver);
462 
463  status = gsl_multimin_test_size (size, this->getTolerance());
464 
465  if(monitor)
466  {
467  gsl_vector* x = gsl_multimin_fminimizer_x(solver);
468  std::vector<double> x_min(dim);
469  for( unsigned int i = 0; i < dim; i++)
470  x_min[i] = gsl_vector_get(x,i);
471 
472  double f = gsl_multimin_fminimizer_minimum(solver);
473 
474  monitor->append( x_min, f, size );
475  }
476 
477  }
478 
479  while ((status == GSL_CONTINUE) && (iter < this->getMaxIterations()));
480 
481  for (unsigned int i = 0; i < dim; i++) {
482  (*m_minimizer)[i] = gsl_vector_get(solver->x, i);
483  }
484 
485  // We're being good human beings and cleaning up the memory we allocated
486  gsl_vector_free(step_size);
487  gsl_multimin_fminimizer_free(solver);
488  gsl_vector_free(x);
489 
490  return;
491 }
GslVector * m_initialPoint
Definition: GslOptimizer.h:157
SolverType m_solver_type
Definition: GslOptimizer.h:160
const BaseScalarFunction< GslVector, GslMatrix > & m_objectiveFunction
Definition: GslOptimizer.h:155
double c_evaluate(const gsl_vector *x, void *context)
Definition: GslOptimizer.C:43
int dim
Definition: ann2fig.cpp:81
double getTolerance() const
Returns the tolerance used to test for an extremum in the optimizer.
Definition: Optimizer.C:57
GslVector m_fstep_size
For use in gradient-free algorithms.
Definition: GslOptimizer.h:163
unsigned int getMaxIterations() const
Returns the maximum number of iterations the optimizer will do.
Definition: Optimizer.C:51
void QUESO::GslOptimizer::minimize_with_gradient ( unsigned int  dim,
OptimizerMonitor monitor 
)
private

Definition at line 293 of file GslOptimizer.C.

References QUESO::OptimizerMonitor::append(), BFGS, BFGS2, QUESO::c_evaluate(), QUESO::c_evaluate_derivative(), QUESO::c_evaluate_with_derivative(), dim, FLETCHER_REEVES_CG, getFdfstepSize(), getLineTolerance(), QUESO::BaseOptimizer::getMaxIterations(), m_initialPoint, m_objectiveFunction, m_solver_type, NELDER_MEAD, NELDER_MEAD2, NELDER_MEAD2_RAND, POLAK_RIBIERE_CG, and STEEPEST_DESCENT.

Referenced by minimize().

294 {
295  // Set initial point
296  gsl_vector * x = gsl_vector_alloc(dim);
297  for (unsigned int i = 0; i < dim; i++) {
298  gsl_vector_set(x, i, (*m_initialPoint)[i]);
299  }
300 
301  // Tell GSL which solver we're using
302  const gsl_multimin_fdfminimizer_type* type = NULL;
303 
304  // We use m_solver_type here because we need the enum
305  switch(m_solver_type)
306  {
307  case(FLETCHER_REEVES_CG):
308  type = gsl_multimin_fdfminimizer_conjugate_fr;
309  break;
310  case(POLAK_RIBIERE_CG):
311  type = gsl_multimin_fdfminimizer_conjugate_pr;
312  break;
313  case(BFGS):
314  type = gsl_multimin_fdfminimizer_vector_bfgs;
315  break;
316  case(BFGS2):
317  type = gsl_multimin_fdfminimizer_vector_bfgs2;
318  break;
319  case(STEEPEST_DESCENT):
320  type = gsl_multimin_fdfminimizer_steepest_descent;
321  break;
322  case(NELDER_MEAD):
323  case(NELDER_MEAD2):
324  case(NELDER_MEAD2_RAND):
325  default:
326  // Wat?!
327  queso_error();
328  }
329 
330  // Init solver
331  gsl_multimin_fdfminimizer * solver =
332  gsl_multimin_fdfminimizer_alloc(type, dim);
333 
334  // Point GSL to the right functions
335  gsl_multimin_function_fdf minusLogPosterior;
336  minusLogPosterior.n = dim;
337  minusLogPosterior.f = &c_evaluate;
338  minusLogPosterior.df = &c_evaluate_derivative;
339  minusLogPosterior.fdf = &c_evaluate_with_derivative;
340  minusLogPosterior.params = (void *)(this);
341 
342  gsl_multimin_fdfminimizer_set(solver, &minusLogPosterior, x, getFdfstepSize(), getLineTolerance());
343 
344  int status;
345  size_t iter = 0;
346 
347  do {
348  iter++;
349  status = gsl_multimin_fdfminimizer_iterate(solver);
350 
351  if (status) {
352  if( m_objectiveFunction.domainSet().env().fullRank() == 0 )
353  {
354  std::cerr << "Error while GSL does optimisation. "
355  << "See below for GSL error type." << std::endl;
356  std::cerr << "Gsl error: " << gsl_strerror(status) << std::endl;
357  }
358  break;
359  }
360 
361  status = gsl_multimin_test_gradient(solver->gradient, this->getTolerance());
362 
363  if(monitor)
364  {
365  gsl_vector* x = gsl_multimin_fdfminimizer_x(solver);
366  std::vector<double> x_min(dim);
367  for( unsigned int i = 0; i < dim; i++)
368  x_min[i] = gsl_vector_get(x,i);
369 
370  double f = gsl_multimin_fdfminimizer_minimum(solver);
371 
372  gsl_vector* grad = gsl_multimin_fdfminimizer_gradient(solver);
373  double grad_norm = gsl_blas_dnrm2(grad);
374 
375  monitor->append( x_min, f, grad_norm );
376  }
377 
378  } while ((status == GSL_CONTINUE) && (iter < this->getMaxIterations()));
379 
380  for (unsigned int i = 0; i < dim; i++) {
381  (*m_minimizer)[i] = gsl_vector_get(solver->x, i);
382  }
383 
384  // We're being good human beings and cleaning up the memory we allocated
385  gsl_multimin_fdfminimizer_free(solver);
386  gsl_vector_free(x);
387 
388  return;
389 }
GslVector * m_initialPoint
Definition: GslOptimizer.h:157
virtual double getLineTolerance() const
Gets the tolerance to use for line minimisation.
Definition: GslOptimizer.C:606
SolverType m_solver_type
Definition: GslOptimizer.h:160
const BaseScalarFunction< GslVector, GslMatrix > & m_objectiveFunction
Definition: GslOptimizer.h:155
void c_evaluate_derivative(const gsl_vector *x, void *context, gsl_vector *derivative)
Definition: GslOptimizer.C:67
void c_evaluate_with_derivative(const gsl_vector *x, void *context, double *f, gsl_vector *derivative)
Definition: GslOptimizer.C:147
double c_evaluate(const gsl_vector *x, void *context)
Definition: GslOptimizer.C:43
int dim
Definition: ann2fig.cpp:81
virtual double getFdfstepSize() const
Gets the step to use in gradient-based solvers.
Definition: GslOptimizer.C:600
unsigned int getMaxIterations() const
Returns the maximum number of iterations the optimizer will do.
Definition: Optimizer.C:51
const GslVector & QUESO::GslOptimizer::minimizer ( ) const

Return the point that minimizes the objective function.

This state will be filled with GSL_NAN if, for some reason, the optimization failed

Definition at line 251 of file GslOptimizer.C.

References m_minimizer.

Referenced by QUESO::StatisticalInverseProblem< P_V, P_M >::solveWithBayesMetropolisHastings().

252 {
253  return *(this->m_minimizer);
254 }
GslVector * m_minimizer
Definition: GslOptimizer.h:158
const BaseScalarFunction< GslVector, GslMatrix > & QUESO::GslOptimizer::objectiveFunction ( ) const

Returns the objective function.

Definition at line 237 of file GslOptimizer.C.

References m_objectiveFunction.

Referenced by QUESO::c_evaluate(), QUESO::c_evaluate_derivative(), and setFstepSize().

238 {
239  return this->m_objectiveFunction;
240 }
const BaseScalarFunction< GslVector, GslMatrix > & m_objectiveFunction
Definition: GslOptimizer.h:155
void QUESO::GslOptimizer::set_line_tol ( double  tol)

Set GSL line minimization tolerance.

Applicable only to gradient-based solvers. Default is 0.1, as recommended by GSL documentation. See GSL documentation for more details.

void QUESO::GslOptimizer::set_solver_type ( SolverType  solver)

Definition at line 256 of file GslOptimizer.C.

References m_solver_type.

Referenced by set_solver_type(), and setSolverType().

257 {
258  queso_deprecated();
259  m_solver_type = solver;
260 }
SolverType m_solver_type
Definition: GslOptimizer.h:160
void QUESO::GslOptimizer::set_solver_type ( std::string &  solver)

Definition at line 549 of file GslOptimizer.C.

References set_solver_type(), and string_to_enum().

550 {
551  queso_deprecated()
552  this->set_solver_type( this->string_to_enum(solver) );
553 }
void set_solver_type(SolverType solver)
Definition: GslOptimizer.C:256
SolverType string_to_enum(std::string &solver)
Definition: GslOptimizer.C:505
void QUESO::GslOptimizer::set_step_size ( const GslVector step_size)

Sets step size used in gradient-free solvers.

By default, the step size used will be a vector of 0.1. Use this method to reset the step_size to the desired values.

Definition at line 493 of file GslOptimizer.C.

References m_fstep_size.

Referenced by setFdfstepSize(), and setFstepSize().

494 {
495  queso_deprecated();
496  m_fstep_size = step_size;
497 }
GslVector m_fstep_size
For use in gradient-free algorithms.
Definition: GslOptimizer.h:163
void QUESO::GslOptimizer::set_step_size ( double  step_size)

Sets step size used in gradient-based solvers.

GSL doesn't document this parameter well, but it seems to be related to the line search, so we default to 1.0 for full step.

Definition at line 499 of file GslOptimizer.C.

References m_fdfstep_size.

500 {
501  queso_deprecated();
502  m_fdfstep_size = step_size;
503 }
double m_fdfstep_size
For use in gradient-based algorithms.
Definition: GslOptimizer.h:166
void QUESO::GslOptimizer::setFdfstepSize ( double  fdfstepSize)
virtual

Sets the step to use in gradient-based solvers.

Reimplemented from QUESO::BaseOptimizer.

Definition at line 575 of file GslOptimizer.C.

References QUESO::BaseOptimizer::m_optionsObj, and set_step_size().

576 {
577  this->m_optionsObj->m_fdfstepSize = fdfstepSize;
578  this->set_step_size(fdfstepSize);
579 }
ScopedPtr< OptimizerOptions >::Type m_optionsObj
Definition: Optimizer.h:126
void set_step_size(const GslVector &step_size)
Sets step size used in gradient-free solvers.
Definition: GslOptimizer.C:493
void QUESO::GslOptimizer::setFstepSize ( double  fstepSize)
virtual

Sets the step size to use in gradient-free solvers.

Reimplemented from QUESO::BaseOptimizer.

Definition at line 563 of file GslOptimizer.C.

References QUESO::GslVector::cwSet(), QUESO::BaseOptimizer::m_optionsObj, objectiveFunction(), and set_step_size().

564 {
565  this->m_optionsObj->m_fstepSize = fstepSize;
566 
567  GslVector fstepSizeVector(
568  objectiveFunction().domainSet().vectorSpace().zeroVector());
569  fstepSizeVector.cwSet(fstepSize);
570 
571  this->set_step_size(fstepSizeVector);
572 }
const BaseScalarFunction< GslVector, GslMatrix > & objectiveFunction() const
Returns the objective function.
Definition: GslOptimizer.C:237
ScopedPtr< OptimizerOptions >::Type m_optionsObj
Definition: Optimizer.h:126
void set_step_size(const GslVector &step_size)
Sets step size used in gradient-free solvers.
Definition: GslOptimizer.C:493
void QUESO::GslOptimizer::setInitialPoint ( const GslVector intialPoint)

Set the point at which the optimization starts.

Definition at line 243 of file GslOptimizer.C.

References m_initialPoint, and QUESO::GslVector::sizeLocal().

Referenced by QUESO::StatisticalInverseProblem< P_V, P_M >::solveWithBayesMetropolisHastings().

244 {
245  for (unsigned int i = 0; i < initialPoint.sizeLocal(); i++) {
246  (*(this->m_initialPoint))[i] = initialPoint[i];
247  }
248 }
GslVector * m_initialPoint
Definition: GslOptimizer.h:157
void QUESO::GslOptimizer::setLineTolerance ( double  lineTolerance)
virtual

Sets the tolerance to use for line minimisation.

Reimplemented from QUESO::BaseOptimizer.

Definition at line 582 of file GslOptimizer.C.

References QUESO::BaseOptimizer::m_optionsObj.

583 {
584  this->m_optionsObj->m_lineTolerance = lineTolerance;
585 }
ScopedPtr< OptimizerOptions >::Type m_optionsObj
Definition: Optimizer.h:126
void QUESO::GslOptimizer::setSolverType ( std::string  solverType)
virtual

Sets the algorithm to use for minimisation.

Reimplemented from QUESO::BaseOptimizer.

Definition at line 556 of file GslOptimizer.C.

References QUESO::BaseOptimizer::m_optionsObj, and set_solver_type().

Referenced by GslOptimizer().

557 {
558  this->m_optionsObj->m_solverType = solverType;
559  this->set_solver_type(solverType);
560 }
ScopedPtr< OptimizerOptions >::Type m_optionsObj
Definition: Optimizer.h:126
void set_solver_type(SolverType solver)
Definition: GslOptimizer.C:256
bool QUESO::GslOptimizer::solver_needs_gradient ( SolverType  solver)
private

Helper function.

Definition at line 262 of file GslOptimizer.C.

References BFGS, BFGS2, FLETCHER_REEVES_CG, NELDER_MEAD, NELDER_MEAD2, NELDER_MEAD2_RAND, POLAK_RIBIERE_CG, and STEEPEST_DESCENT.

Referenced by minimize().

263 {
264  bool gradient_needed = false;
265 
266  switch(solver)
267  {
268  case(FLETCHER_REEVES_CG):
269  case(POLAK_RIBIERE_CG):
270  case(BFGS):
271  case(BFGS2):
272  case(STEEPEST_DESCENT):
273  {
274  gradient_needed = true;
275  break;
276  }
277  case(NELDER_MEAD):
278  case(NELDER_MEAD2):
279  case(NELDER_MEAD2_RAND):
280  {
281  break;
282  }
283  default:
284  {
285  // Wat?!
286  queso_error();
287  }
288  } // switch(solver)
289 
290  return gradient_needed;
291 }
GslOptimizer::SolverType QUESO::GslOptimizer::string_to_enum ( std::string &  solver)

Definition at line 505 of file GslOptimizer.C.

References BFGS, BFGS2, FLETCHER_REEVES_CG, m_objectiveFunction, NELDER_MEAD, NELDER_MEAD2, NELDER_MEAD2_RAND, POLAK_RIBIERE_CG, and STEEPEST_DESCENT.

Referenced by set_solver_type().

506 {
507  SolverType solver_type;
508 
509  if( solver == std::string("fletcher_reeves_cg") )
510  solver_type = FLETCHER_REEVES_CG;
511  else if( solver == std::string("polak_ribiere_cg") )
512  solver_type = POLAK_RIBIERE_CG;
513  else if( solver == std::string("bfgs") )
514  solver_type = BFGS;
515  else if( solver == std::string("bfgs2") )
516  solver_type = BFGS2;
517  else if( solver == std::string("steepest_decent") ) {
518  queso_deprecated();
519  solver_type = STEEPEST_DESCENT;
520  }
521  else if( solver == std::string("steepest_descent") )
522  solver_type = STEEPEST_DESCENT;
523  else if( solver == std::string("nelder_mead") )
524  solver_type = NELDER_MEAD;
525  else if( solver == std::string("nelder_mead2") )
526  solver_type = NELDER_MEAD2;
527  else if( solver == std::string("nelder_mead2_rand") )
528  solver_type = NELDER_MEAD2_RAND;
529  else
530  {
531  if( m_objectiveFunction.domainSet().env().fullRank() == 0 )
532  {
533  std::cerr << "Error: Invalid GslOptimizer solver name: " << solver << std::endl
534  << " Valids choices are: fletcher_reeves_cg" << std::endl
535  << " polak_ribiere_cg" << std::endl
536  << " bfgs" << std::endl
537  << " bfgs2" << std::endl
538  << " steepest_descent" << std::endl
539  << " nelder_mead" << std::endl
540  << " nelder_mead2" << std::endl
541  << " nelder_mead2_rand" << std::endl;
542  }
543  queso_error();
544  }
545 
546  return solver_type;
547 }
const BaseScalarFunction< GslVector, GslMatrix > & m_objectiveFunction
Definition: GslOptimizer.h:155

Member Data Documentation

double QUESO::GslOptimizer::m_fdfstep_size
private

For use in gradient-based algorithms.

Definition at line 166 of file GslOptimizer.h.

Referenced by set_step_size().

GslVector QUESO::GslOptimizer::m_fstep_size
private

For use in gradient-free algorithms.

Definition at line 163 of file GslOptimizer.h.

Referenced by GslOptimizer(), minimize_no_gradient(), and set_step_size().

GslVector* QUESO::GslOptimizer::m_initialPoint
private
double QUESO::GslOptimizer::m_line_tol
private

Line minimization tolerance in gradient-based algorithms.

Definition at line 169 of file GslOptimizer.h.

GslVector* QUESO::GslOptimizer::m_minimizer
private

Definition at line 158 of file GslOptimizer.h.

Referenced by GslOptimizer(), and minimizer().

const BaseScalarFunction<GslVector, GslMatrix>& QUESO::GslOptimizer::m_objectiveFunction
private
SolverType QUESO::GslOptimizer::m_solver_type
private

The documentation for this class was generated from the following files:

Generated on Sat Apr 22 2017 14:04:39 for queso-0.57.0 by  doxygen 1.8.5