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