queso-0.53.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...
 
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...
 
- Public Member Functions inherited from QUESO::BaseOptimizer
 BaseOptimizer ()
 Default constructor. 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
 

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 157 of file GslOptimizer.C.

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

159  : BaseOptimizer(),
161  m_initialPoint(new GslVector(objectiveFunction.domainSet().
162  vectorSpace().zeroVector())),
163  m_minimizer(new GslVector(this->m_objectiveFunction.domainSet().
164  vectorSpace().zeroVector())),
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 }
GslVector * m_minimizer
Definition: GslOptimizer.h:123
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
BaseOptimizer()
Default constructor.
Definition: Optimizer.C:29
GslVector * m_initialPoint
Definition: GslOptimizer.h:122
void cwSet(double value)
Component-wise sets all values to this with value.
Definition: GslVector.C:339
const BaseScalarFunction< GslVector, GslMatrix > & m_objectiveFunction
Definition: GslOptimizer.h:120
const BaseScalarFunction< GslVector, GslMatrix > & objectiveFunction() const
Returns the objective function.
Definition: GslOptimizer.C:211
SolverType m_solver_type
Definition: GslOptimizer.h:125
double m_fdfstep_size
For use in gradient-based algorithms.
Definition: GslOptimizer.h:131
QUESO::GslOptimizer::~GslOptimizer ( )
virtual

Destructor.

Definition at line 177 of file GslOptimizer.C.

References m_initialPoint.

178 {
179  delete this->m_initialPoint;
180 }
GslVector * m_initialPoint
Definition: GslOptimizer.h:122

Member Function Documentation

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 182 of file GslOptimizer.C.

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

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

182  {
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 }
GslVector * m_initialPoint
Definition: GslOptimizer.h:122
#define queso_error()
Definition: asserts.h:53
bool solver_needs_gradient(SolverType solver)
Helper function.
Definition: GslOptimizer.C:235
const BaseScalarFunction< GslVector, GslMatrix > & m_objectiveFunction
Definition: GslOptimizer.h:120
void minimize_no_gradient(unsigned int dim, OptimizerMonitor *monitor)
Definition: GslOptimizer.C:363
SolverType m_solver_type
Definition: GslOptimizer.h:125
void minimize_with_gradient(unsigned int dim, OptimizerMonitor *monitor)
Definition: GslOptimizer.C:266
int dim
Definition: ann2fig.cpp:81
void QUESO::GslOptimizer::minimize_no_gradient ( unsigned int  dim,
OptimizerMonitor monitor 
)
private

Definition at line 363 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_error, and STEEPEST_DESCENT.

Referenced by minimize().

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  }
GslVector m_fstep_size
For use in gradient-free algorithms.
Definition: GslOptimizer.h:128
GslVector * m_initialPoint
Definition: GslOptimizer.h:122
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
double getTolerance() const
Returns the tolerance used to test for an extremum in the optimizer.
Definition: Optimizer.C:47
const BaseScalarFunction< GslVector, GslMatrix > & m_objectiveFunction
Definition: GslOptimizer.h:120
SolverType m_solver_type
Definition: GslOptimizer.h:125
int dim
Definition: ann2fig.cpp:81
void QUESO::GslOptimizer::minimize_with_gradient ( unsigned int  dim,
OptimizerMonitor monitor 
)
private

Definition at line 266 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, QUESO::BaseOptimizer::getMaxIterations(), m_fdfstep_size, m_initialPoint, m_line_tol, m_objectiveFunction, m_solver_type, NELDER_MEAD, NELDER_MEAD2, NELDER_MEAD2_RAND, POLAK_RIBIERE_CG, queso_error, and STEEPEST_DESCENT.

Referenced by minimize().

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  }
double m_line_tol
Line minimization tolerance in gradient-based algorithms.
Definition: GslOptimizer.h:134
void c_evaluate_derivative(const gsl_vector *x, void *context, gsl_vector *derivative)
Definition: GslOptimizer.C:67
GslVector * m_initialPoint
Definition: GslOptimizer.h:122
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
void c_evaluate_with_derivative(const gsl_vector *x, void *context, double *f, gsl_vector *derivative)
Definition: GslOptimizer.C:149
const BaseScalarFunction< GslVector, GslMatrix > & m_objectiveFunction
Definition: GslOptimizer.h:120
SolverType m_solver_type
Definition: GslOptimizer.h:125
double m_fdfstep_size
For use in gradient-based algorithms.
Definition: GslOptimizer.h:131
int dim
Definition: ann2fig.cpp:81
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 225 of file GslOptimizer.C.

References m_minimizer.

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

226 {
227  return *(this->m_minimizer);
228 }
GslVector * m_minimizer
Definition: GslOptimizer.h:123
const BaseScalarFunction< GslVector, GslMatrix > & QUESO::GslOptimizer::objectiveFunction ( ) const

Returns the objective function.

Definition at line 211 of file GslOptimizer.C.

References m_objectiveFunction.

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

212 {
213  return this->m_objectiveFunction;
214 }
const BaseScalarFunction< GslVector, GslMatrix > & m_objectiveFunction
Definition: GslOptimizer.h:120
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 230 of file GslOptimizer.C.

References m_solver_type.

Referenced by set_solver_type().

231  {
232  m_solver_type = solver;
233  }
SolverType m_solver_type
Definition: GslOptimizer.h:125
void QUESO::GslOptimizer::set_solver_type ( std::string &  solver)

Definition at line 515 of file GslOptimizer.C.

References set_solver_type(), and string_to_enum().

516  {
517  this->set_solver_type( this->string_to_enum(solver) );
518  }
void set_solver_type(SolverType solver)
Definition: GslOptimizer.C:230
SolverType string_to_enum(std::string &solver)
Definition: GslOptimizer.C:475
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 465 of file GslOptimizer.C.

References m_fstep_size.

466  {
467  m_fstep_size = step_size;
468  }
GslVector m_fstep_size
For use in gradient-free algorithms.
Definition: GslOptimizer.h:128
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 470 of file GslOptimizer.C.

References m_fdfstep_size.

471  {
472  m_fdfstep_size = step_size;
473  }
double m_fdfstep_size
For use in gradient-based algorithms.
Definition: GslOptimizer.h:131
void QUESO::GslOptimizer::setInitialPoint ( const GslVector intialPoint)

Set the point at which the optimization starts.

Definition at line 217 of file GslOptimizer.C.

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

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

218 {
219  for (unsigned int i = 0; i < initialPoint.sizeLocal(); i++) {
220  (*(this->m_initialPoint))[i] = initialPoint[i];
221  }
222 }
GslVector * m_initialPoint
Definition: GslOptimizer.h:122
bool QUESO::GslOptimizer::solver_needs_gradient ( SolverType  solver)
private

Helper function.

Definition at line 235 of file GslOptimizer.C.

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

Referenced by minimize().

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  }
#define queso_error()
Definition: asserts.h:53
GslOptimizer::SolverType QUESO::GslOptimizer::string_to_enum ( std::string &  solver)

Definition at line 475 of file GslOptimizer.C.

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

Referenced by set_solver_type().

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  }
#define queso_error()
Definition: asserts.h:53
const BaseScalarFunction< GslVector, GslMatrix > & m_objectiveFunction
Definition: GslOptimizer.h:120

Member Data Documentation

double QUESO::GslOptimizer::m_fdfstep_size
private

For use in gradient-based algorithms.

Definition at line 131 of file GslOptimizer.h.

Referenced by minimize_with_gradient(), and set_step_size().

GslVector QUESO::GslOptimizer::m_fstep_size
private

For use in gradient-free algorithms.

Definition at line 128 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 134 of file GslOptimizer.h.

Referenced by minimize_with_gradient().

GslVector* QUESO::GslOptimizer::m_minimizer
private

Definition at line 123 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 Thu Jun 11 2015 13:52:34 for queso-0.53.0 by  doxygen 1.8.5