bpp-core  2.2.0
ConjugateGradientMultiDimensions.cpp
Go to the documentation of this file.
1 //
2 // File: ConjugateGradientMultiDimensions.cpp
3 // Created by: Julien Dutheil
4 // Created on: Wed Apr 11 16:51 2007
5 //
6 
7 /*
8 Copyright or © or Copr. Bio++ Development Team, (November 19, 2004)
9 
10 This software is a computer program whose purpose is to provide classes
11 for numerical calculus.
12 
13 This software is governed by the CeCILL license under French law and
14 abiding by the rules of distribution of free software. You can use,
15 modify and/ or redistribute the software under the terms of the CeCILL
16 license as circulated by CEA, CNRS and INRIA at the following URL
17 "http://www.cecill.info".
18 
19 As a counterpart to the access to the source code and rights to copy,
20 modify and redistribute granted by the license, users are provided only
21 with a limited warranty and the software's author, the holder of the
22 economic rights, and the successive licensors have only limited
23 liability.
24 
25 In this respect, the user's attention is drawn to the risks associated
26 with loading, using, modifying and/or developing or reproducing the
27 software by the user in light of its specific status of free software,
28 that may mean that it is complicated to manipulate, and that also
29 therefore means that it is reserved for developers and experienced
30 professionals having in-depth computer knowledge. Users are therefore
31 encouraged to load and test the software's suitability as regards their
32 requirements in conditions enabling the security of their systems and/or
33 data to be ensured and, more generally, to use and operate it in the
34 same conditions as regards security.
35 
36 The fact that you are presently reading this means that you have had
37 knowledge of the CeCILL license and that you accept its terms.
38 */
39 
42 
43 using namespace bpp;
44 
45 /******************************************************************************/
46 
48  AbstractOptimizer(function), optimizer_(function),
49  xi_(), h_(), g_(), f1dim_(function)
50 {
53 }
54 
55 /******************************************************************************/
56 
58 {
59  size_t nbParams = params.size();
60  g_.resize(nbParams);
61  h_.resize(nbParams);
62  xi_.resize(nbParams);
63  getFunction_()->enableFirstOrderDerivatives(true);
64  getFunction_()->setParameters(params);
65  getGradient(xi_);
66  for(size_t i = 0; i < nbParams; i++)
67  {
68  g_[i] = -xi_[i];
69  xi_[i] = h_[i] = g_[i];
70  }
71 }
72 
73 /******************************************************************************/
74 
76 {
77  double gg, gam, f, dgg;
78  size_t n = getParameters().size();
79  //Loop over iterations.
82  getParameters_(), xi_, getStopCondition()->getTolerance(),
83  0, 0, getVerbose() > 0 ? getVerbose() - 1 : 0);
85  f = getFunction()->f(getParameters());
86  if (tolIsReached_)
87  {
88  return f;
89  }
91  dgg = gg = 0.0;
92  for (unsigned j = 0; j < n; j++)
93  {
94  gg += g_[j] * g_[j];
95  /* dgg += xi[j] * xi[j]; */ //This statement for Fletcher-Reeves.
96  dgg += (xi_[j] + g_[j]) * xi_[j]; //This statement for Polak-Ribiere.
97  }
98  if (gg == 0.0)
99  {
100  //Unlikely. If gradient is exactly zero then
101  return f;
102  }
103  gam = dgg / gg;
104 
105  if (!(std::isnan(gam) || std::isinf(gam)))
106  {
107  for(unsigned int j = 0; j < n; j++)
108  {
109  g_[j] = -xi_[j];
110  xi_[j] = h_[j] = g_[j] + gam * h_[j];
111  }
112  }
113 
114  return f;
115 }
116 
117 /******************************************************************************/
118 
119 void ConjugateGradientMultiDimensions::getGradient(std::vector<double>& gradient) const
120 {
121  for (size_t i = 0; i < gradient.size(); ++i)
122  {
123  gradient[i] = getFunction()->getFirstOrderDerivative(getParameters()[i].getName());
124  }
125 }
126 
127 /******************************************************************************/
128 
virtual double f(const ParameterList &parameters)
Get the value of the function according to a given set of parameters.
Definition: Functions.h:117
OptimizationStopCondition * getStopCondition()
Get the stop condition of the optimization algorithm.
unsigned int getVerbose() const
Get the verbose level.
This class allows to perform a correspondence analysis.
size_t size() const
Definition: ParameterList.h:90
ParameterList & getParameters_()
double doStep()
This function is called by the step() method and contains all calculations.
Stop condition on function value.
const DerivableFirstOrder * getFunction() const
Get the current function being optimized.
virtual double getFirstOrderDerivative(const std::string &variable) const =0
Get the derivative of the function at the current point.
The parameter list object.
Definition: ParameterList.h:61
virtual void enableFirstOrderDerivatives(bool yn)=0
Tell if derivatives must be computed.
bool tolIsReached_
Tell if the tolerance level has been reached.
void setDefaultStopCondition_(OptimizationStopCondition *osc)
const ParameterList & getParameters() const
ConjugateGradientMultiDimensions(DerivableFirstOrder *function)
void doInit(const ParameterList &params)
This function is called by the init() method and contains all calculations.
Exception base class.
Definition: Exceptions.h:57
unsigned int nbEval_
The current number of function evaluations achieved.
Partial implementation of the Optimizer interface.
void setStopCondition(const OptimizationStopCondition &stopCondition)
Set the stop condition of the optimization algorithm.
This is the abstract class for first order derivable functions.
Definition: Functions.h:129
OptimizationStopCondition * getDefaultStopCondition()
Get the default stop condition of the optimization algorithm.
static unsigned int lineMinimization(DirectionFunction &f1dim, ParameterList &parameters, std::vector< double > &xi, double tolerance, OutputStream *profiler=0, OutputStream *messenger=0, unsigned int verbose=2)
void getGradient(std::vector< double > &gradient) const