bpp-phyl  2.2.0
OptimizationTools.h
Go to the documentation of this file.
1 //
2 // File: OptimizationTools.h
3 // Created by: Julien Dutheil
4 // Created on: Sun Dec 14 09:43:32 2003
5 //
6 
7 /*
8  Copyright or © or Copr. Bio++ Development Team, (November 16, 2004)
9 
10  This software is a computer program whose purpose is to provide classes
11  for phylogenetic data analysis.
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 
40 #ifndef _OPTIMIZATIONTOOLS_H_
41 #define _OPTIMIZATIONTOOLS_H_
42 
46 #include "NNITopologySearch.h"
48 #include "TreeTemplate.h"
51 
52 #include <Bpp/Io/OutputStream.h>
53 #include <Bpp/App/ApplicationTools.h>
54 #include <Bpp/Numeric/Function/SimpleNewtonMultiDimensions.h>
55 
56 namespace bpp
57 {
58 
62 class NaNListener: public OptimizationListener
63 {
64  private:
65  Optimizer* optimizer_;
66  Function* function_;
67 
68  public:
69  NaNListener(Optimizer* optimizer, Function* function): optimizer_(optimizer), function_(function) {}
70 
74  {}
75 
77  {
79  function_ = lr.function_;
80  return *this;
81  }
82 
83  public:
84  void optimizationInitializationPerformed(const OptimizationEvent &event) {}
85  void optimizationStepPerformed(const OptimizationEvent &event) throw (Exception)
86  {
87  if (std::isnan(optimizer_->getFunction()->getValue()))
88  {
89  cerr << "Oups... something abnormal happened!" << endl;
90  function_->getParameters().printParameters(cerr);
91  throw Exception("Optimization failed because likelihood function returned NaN.");
92  }
93  }
94  bool listenerModifiesParameters () const { return false; }
95 };
96 
97 
98 
103  public virtual TopologyListener
104 {
105 private:
107  ParameterList parameters_;
108  double tolerance_;
109  OutputStream* messenger_;
110  OutputStream* profiler_;
111  unsigned int verbose_;
112  unsigned int optimizeCounter_;
113  unsigned int optimizeNumerical_;
114  std::string optMethod_;
115  unsigned int nStep_;
117 
118 public:
137  NNITopologySearch* ts,
138  const ParameterList& parameters,
139  double tolerance,
140  OutputStream* messenger,
141  OutputStream* profiler,
142  unsigned int verbose,
143  const std::string& optMethod,
144  unsigned int nStep,
145  bool reparametrization) :
146  topoSearch_(ts),
147  parameters_(parameters),
148  tolerance_(tolerance),
149  messenger_(messenger),
150  profiler_(profiler),
151  verbose_(verbose),
152  optimizeCounter_(0),
154  optMethod_(optMethod),
155  nStep_(nStep),
156  reparametrization_(reparametrization) {}
157 
163  profiler_(tl.profiler_),
164  verbose_(tl.verbose_),
168  nStep_(tl.nStep_),
170  {}
171 
173  {
176  tolerance_ = tl.tolerance_;
177  messenger_ = tl.messenger_;
178  profiler_ = tl.profiler_;
179  verbose_ = tl.verbose_;
182  optMethod_ = tl.optMethod_;
183  nStep_ = tl.nStep_;
185  return *this;
186  }
187 
188  NNITopologyListener* clone() const { return new NNITopologyListener(*this); }
189 
190  virtual ~NNITopologyListener() {}
191 
192 public:
196 };
197 
202  public TopologyListener
203 {
204 private:
206  ParameterList parameters_;
207  double tolerance_;
208  OutputStream* messenger_;
209  OutputStream* profiler_;
210  unsigned int verbose_;
211  unsigned int optimizeCounter_;
212  unsigned int optimizeNumerical_;
213  std::string optMethod_;
215 
216 public:
234  NNITopologySearch* ts,
235  const ParameterList& parameters,
236  double tolerance,
237  OutputStream* messenger,
238  OutputStream* profiler,
239  unsigned int verbose,
240  const std::string& optMethod,
241  bool reparametrization) :
242  topoSearch_(ts),
243  parameters_(parameters),
244  tolerance_(tolerance),
245  messenger_(messenger),
246  profiler_(profiler),
247  verbose_(verbose),
248  optimizeCounter_(0),
250  optMethod_(optMethod),
251  reparametrization_(reparametrization) {}
252 
258  profiler_(tl.profiler_),
259  verbose_(tl.verbose_),
264  {}
265 
267  {
270  tolerance_ = tl.tolerance_;
271  messenger_ = tl.messenger_;
272  profiler_ = tl.profiler_;
273  verbose_ = tl.verbose_;
276  optMethod_ = tl.optMethod_;
278  return *this;
279  }
280 
281  NNITopologyListener2* clone() const { return new NNITopologyListener2(*this); }
282 
284 
285 public:
289 };
290 
291 
301 {
302 public:
304  virtual ~OptimizationTools();
305 
306 public:
307  static std::string OPTIMIZATION_GRADIENT;
308  static std::string OPTIMIZATION_NEWTON;
309  static std::string OPTIMIZATION_BRENT;
310  static std::string OPTIMIZATION_BFGS;
311 
338  static unsigned int optimizeNumericalParameters(
340  const ParameterList& parameters,
341  OptimizationListener* listener = 0,
342  unsigned int nstep = 1,
343  double tolerance = 0.000001,
344  unsigned int tlEvalMax = 1000000,
345  OutputStream* messageHandler = ApplicationTools::message,
346  OutputStream* profiler = ApplicationTools::message,
347  bool reparametrization = false,
348  unsigned int verbose = 1,
349  const std::string& optMethodDeriv = OPTIMIZATION_NEWTON,
350  const std::string& optMethodModel = OPTIMIZATION_BRENT)
351  throw (Exception);
352 
375  static unsigned int optimizeNumericalParameters2(
377  const ParameterList& parameters,
378  OptimizationListener* listener = 0,
379  double tolerance = 0.000001,
380  unsigned int tlEvalMax = 1000000,
381  OutputStream* messageHandler = ApplicationTools::message,
382  OutputStream* profiler = ApplicationTools::message,
383  bool reparametrization = false,
384  bool useClock = false,
385  unsigned int verbose = 1,
386  const std::string& optMethodDeriv = OPTIMIZATION_NEWTON)
387  throw (Exception);
388 
410  static unsigned int optimizeBranchLengthsParameters(
412  const ParameterList& parameters,
413  OptimizationListener* listener = 0,
414  double tolerance = 0.000001,
415  unsigned int tlEvalMax = 1000000,
416  OutputStream* messageHandler = ApplicationTools::message,
417  OutputStream* profiler = ApplicationTools::message,
418  unsigned int verbose = 1,
419  const std::string& optMethodDeriv = OPTIMIZATION_NEWTON)
420  throw (Exception);
421 
447  static unsigned int optimizeNumericalParametersWithGlobalClock(
449  const ParameterList& parameters,
450  OptimizationListener* listener = 0,
451  unsigned int nstep = 1,
452  double tolerance = 0.000001,
453  unsigned int tlEvalMax = 1000000,
454  OutputStream* messageHandler = ApplicationTools::message,
455  OutputStream* profiler = ApplicationTools::message,
456  unsigned int verbose = 1,
457  const std::string& optMethodDeriv = OPTIMIZATION_GRADIENT)
458  throw (Exception);
459 
482  const ParameterList& parameters,
483  OptimizationListener* listener = 0,
484  double tolerance = 0.000001,
485  unsigned int tlEvalMax = 1000000,
486  OutputStream* messageHandler = ApplicationTools::message,
487  OutputStream* profiler = ApplicationTools::message,
488  unsigned int verbose = 1,
489  const std::string& optMethodDeriv = OPTIMIZATION_GRADIENT)
490  throw (Exception);
491 
492 private:
494  public virtual Function,
495  public ParametrizableAdapter
496  {
497 private:
499  mutable ParameterList brLen_, lambda_;
500 
501 public:
503 
505  tl_(sf.tl_),
506  brLen_(sf.brLen_),
507  lambda_(sf.lambda_)
508  {}
509 
511  {
512  tl_ = sf.tl_;
513  brLen_ = sf.brLen_;
514  lambda_ = sf.lambda_;
515  return *this;
516  }
517 
518  virtual ~ScaleFunction();
519 
520 #ifndef NO_VIRTUAL_COV
522 #else
523  Clonable*
524 #endif
525  clone() const { return new ScaleFunction(*this); }
526 
527 public:
528  void setParameters(const ParameterList& lambda) throw (ParameterNotFoundException, ConstraintException);
529  double getValue() const throw (ParameterException);
530  const ParameterList& getParameters() const throw (Exception) { return lambda_; }
531  const Parameter& getParameter(const std::string& name) const throw (ParameterNotFoundException)
532  {
533  if (name == "lambda") return lambda_[0];
534  else throw ParameterNotFoundException("ScaleFunction::getParameter.", name);
535  }
536  double getParameterValue(const std::string& name) const throw (ParameterNotFoundException)
537  {
538  return lambda_.getParameter(name).getValue();
539  }
540  size_t getNumberOfParameters() const { return 1; }
541  size_t getNumberOfIndependentParameters() const { return 1; }
542  };
543 
544 public:
567  static unsigned int optimizeTreeScale(
568  TreeLikelihood* tl,
569  double tolerance = 0.000001,
570  unsigned int tlEvalMax = 1000000,
571  OutputStream* messageHandler = ApplicationTools::message,
572  OutputStream* profiler = ApplicationTools::message,
573  unsigned int verbose = 1)
574  throw (Exception);
575 
617  const ParameterList& parameters,
618  bool optimizeNumFirst = true,
619  double tolBefore = 100,
620  double tolDuring = 100,
621  unsigned int tlEvalMax = 1000000,
622  unsigned int numStep = 1,
623  OutputStream* messageHandler = ApplicationTools::message,
624  OutputStream* profiler = ApplicationTools::message,
625  bool reparametrization = false,
626  unsigned int verbose = 1,
627  const std::string& optMethod = OptimizationTools::OPTIMIZATION_NEWTON,
628  unsigned int nStep = 1,
629  const std::string& nniMethod = NNITopologySearch::PHYML)
630  throw (Exception);
631 
672  const ParameterList& parameters,
673  bool optimizeNumFirst = true,
674  double tolBefore = 100,
675  double tolDuring = 100,
676  unsigned int tlEvalMax = 1000000,
677  unsigned int numStep = 1,
678  OutputStream* messageHandler = ApplicationTools::message,
679  OutputStream* profiler = ApplicationTools::message,
680  bool reparametrization = false,
681  unsigned int verbose = 1,
682  const std::string& optMethod = OptimizationTools::OPTIMIZATION_NEWTON,
683  const std::string& nniMethod = NNITopologySearch::PHYML)
684  throw (Exception);
685 
702  unsigned int verbose = 1);
703 
720  static DistanceMatrix* estimateDistanceMatrix(
721  DistanceEstimation& estimationMethod,
722  const ParameterList& parametersToIgnore,
723  const std::string& param = DISTANCEMETHOD_INIT,
724  unsigned int verbose = 0) throw (Exception);
725 
751  DistanceEstimation& estimationMethod,
752  AgglomerativeDistanceMethod& reconstructionMethod,
753  const ParameterList& parametersToIgnore,
754  bool optimizeBrLen = false,
755  const std::string& param = DISTANCEMETHOD_INIT,
756  double tolerance = 0.000001,
757  unsigned int tlEvalMax = 1000000,
758  OutputStream* profiler = 0,
759  OutputStream* messenger = 0,
760  unsigned int verbose = 0) throw (Exception);
761 
762 public:
763  static std::string DISTANCEMETHOD_INIT;
764  static std::string DISTANCEMETHOD_PAIRWISE;
766 };
767 } // end of namespace bpp.
768 
769 #endif // _OPTIMIZATIONTOOLS_H_
770 
static NNIHomogeneousTreeLikelihood * optimizeTreeNNI(NNIHomogeneousTreeLikelihood *tl, const ParameterList &parameters, bool optimizeNumFirst=true, double tolBefore=100, double tolDuring=100, unsigned int tlEvalMax=1000000, unsigned int numStep=1, OutputStream *messageHandler=ApplicationTools::message, OutputStream *profiler=ApplicationTools::message, bool reparametrization=false, unsigned int verbose=1, const std::string &optMethod=OptimizationTools::OPTIMIZATION_NEWTON, unsigned int nStep=1, const std::string &nniMethod=NNITopologySearch::PHYML)
Optimize all parameters from a TreeLikelihood object, including tree topology using Nearest Neighbor ...
void setNumericalOptimizationCounter(unsigned int c)
static unsigned int optimizeBranchLengthsParameters(DiscreteRatesAcrossSitesTreeLikelihood *tl, const ParameterList &parameters, OptimizationListener *listener=0, double tolerance=0.000001, unsigned int tlEvalMax=1000000, OutputStream *messageHandler=ApplicationTools::message, OutputStream *profiler=ApplicationTools::message, unsigned int verbose=1, const std::string &optMethodDeriv=OPTIMIZATION_NEWTON)
Optimize branch lengths parameters of a TreeLikelihood function.
static unsigned int optimizeNumericalParametersWithGlobalClock2(DiscreteRatesAcrossSitesClockTreeLikelihood *cl, const ParameterList &parameters, OptimizationListener *listener=0, double tolerance=0.000001, unsigned int tlEvalMax=1000000, OutputStream *messageHandler=ApplicationTools::message, OutputStream *profiler=ApplicationTools::message, unsigned int verbose=1, const std::string &optMethodDeriv=OPTIMIZATION_GRADIENT)
Optimize numerical parameters assuming a global clock (branch heights, substitution model & rate dist...
static DistanceMatrix * estimateDistanceMatrix(DistanceEstimation &estimationMethod, const ParameterList &parametersToIgnore, const std::string &param=DISTANCEMETHOD_INIT, unsigned int verbose=0)
Estimate a distance matrix using maximum likelihood.
void topologyChangeTested(const TopologyChangeEvent &event)
Notify a topology change event.
NNITopologyListener * clone() const
static std::string DISTANCEMETHOD_PAIRWISE
Listener used internally by the optimizeTreeNNI2 method.
NNITopologyListener2(NNITopologySearch *ts, const ParameterList &parameters, double tolerance, OutputStream *messenger, OutputStream *profiler, unsigned int verbose, const std::string &optMethod, bool reparametrization)
Build a new NNITopologyListener2 object.
NNITopologyListener2(const NNITopologyListener2 &tl)
NaNListener(const NaNListener &lr)
void setParameters(const ParameterList &lambda)
The TreeLikelihood interface.
A listener which capture NaN function values and throw an exception in case this happens.
Double recursive implementation of interface TreeParsimonyScore.
static std::string DISTANCEMETHOD_INIT
STL namespace.
The phylogenetic tree class.
static TreeTemplate< Node > * buildDistanceTree(DistanceEstimation &estimationMethod, AgglomerativeDistanceMethod &reconstructionMethod, const ParameterList &parametersToIgnore, bool optimizeBrLen=false, const std::string &param=DISTANCEMETHOD_INIT, double tolerance=0.000001, unsigned int tlEvalMax=1000000, OutputStream *profiler=0, OutputStream *messenger=0, unsigned int verbose=0)
Build a tree using a distance method.
static std::string OPTIMIZATION_GRADIENT
This class adds support for NNI topology estimation to the DRHomogeneousTreeLikelihood class...
void topologyChangeTested(const TopologyChangeEvent &event)
Notify a topology change event.
static std::string DISTANCEMETHOD_ITERATIONS
Optimizer * optimizer_
NaNListener & operator=(const NaNListener &lr)
void optimizationStepPerformed(const OptimizationEvent &event)
NNITopologyListener & operator=(const NNITopologyListener &tl)
Interface for agglomerative distance methods.
const Parameter & getParameter(const std::string &name) const
NNITopologyListener(const NNITopologyListener &tl)
NNITopologyListener2 * clone() const
ScaleFunction & operator=(const ScaleFunction &sf)
void optimizationInitializationPerformed(const OptimizationEvent &event)
static unsigned int optimizeNumericalParameters(DiscreteRatesAcrossSitesTreeLikelihood *tl, const ParameterList &parameters, OptimizationListener *listener=0, unsigned int nstep=1, double tolerance=0.000001, unsigned int tlEvalMax=1000000, OutputStream *messageHandler=ApplicationTools::message, OutputStream *profiler=ApplicationTools::message, bool reparametrization=false, unsigned int verbose=1, const std::string &optMethodDeriv=OPTIMIZATION_NEWTON, const std::string &optMethodModel=OPTIMIZATION_BRENT)
Optimize numerical parameters (branch length, substitution model & rate distribution) of a TreeLikeli...
Interface for rate across sites (RAS) implementation.
static std::string OPTIMIZATION_BRENT
static unsigned int optimizeNumericalParameters2(DiscreteRatesAcrossSitesTreeLikelihood *tl, const ParameterList &parameters, OptimizationListener *listener=0, double tolerance=0.000001, unsigned int tlEvalMax=1000000, OutputStream *messageHandler=ApplicationTools::message, OutputStream *profiler=ApplicationTools::message, bool reparametrization=false, bool useClock=false, unsigned int verbose=1, const std::string &optMethodDeriv=OPTIMIZATION_NEWTON)
Optimize numerical parameters (branch length, substitution model & rate distribution) of a TreeLikeli...
void setNumericalOptimizationCounter(unsigned int c)
Interface for likelihood computation with a global clock and rate across sites variation.
Optimization methods for phylogenetic inference.
NNITopologyListener(NNITopologySearch *ts, const ParameterList &parameters, double tolerance, OutputStream *messenger, OutputStream *profiler, unsigned int verbose, const std::string &optMethod, unsigned int nStep, bool reparametrization)
Build a new NNITopologyListener object.
Class for notifying new toplogy change events.
NNITopologySearch * topoSearch_
void topologyChangeSuccessful(const TopologyChangeEvent &event)
Tell that a topology change is definitive.
The phylogenetic node class.
Definition: Node.h:90
NNITopologySearch * topoSearch_
const ParameterList & getParameters() const
bool listenerModifiesParameters() const
static NNIHomogeneousTreeLikelihood * optimizeTreeNNI2(NNIHomogeneousTreeLikelihood *tl, const ParameterList &parameters, bool optimizeNumFirst=true, double tolBefore=100, double tolDuring=100, unsigned int tlEvalMax=1000000, unsigned int numStep=1, OutputStream *messageHandler=ApplicationTools::message, OutputStream *profiler=ApplicationTools::message, bool reparametrization=false, unsigned int verbose=1, const std::string &optMethod=OptimizationTools::OPTIMIZATION_NEWTON, const std::string &nniMethod=NNITopologySearch::PHYML)
Optimize all parameters from a TreeLikelihood object, including tree topology using Nearest Neighbor ...
ScaleFunction(const ScaleFunction &sf)
NNITopologyListener2 & operator=(const NNITopologyListener2 &tl)
void topologyChangeSuccessful(const TopologyChangeEvent &event)
Tell that a topology change is definitive.
Estimate a distance matrix from sequence data, according to a given model.
static std::string OPTIMIZATION_BFGS
Implement this interface to be notified when the topology of a tree has changed during topology searc...
static unsigned int optimizeTreeScale(TreeLikelihood *tl, double tolerance=0.000001, unsigned int tlEvalMax=1000000, OutputStream *messageHandler=ApplicationTools::message, OutputStream *profiler=ApplicationTools::message, unsigned int verbose=1)
Optimize the scale of a TreeLikelihood.
Listener used internally by the optimizeTreeNNI method.
static std::string OPTIMIZATION_NEWTON
double getParameterValue(const std::string &name) const
static unsigned int optimizeNumericalParametersWithGlobalClock(DiscreteRatesAcrossSitesClockTreeLikelihood *cl, const ParameterList &parameters, OptimizationListener *listener=0, unsigned int nstep=1, double tolerance=0.000001, unsigned int tlEvalMax=1000000, OutputStream *messageHandler=ApplicationTools::message, OutputStream *profiler=ApplicationTools::message, unsigned int verbose=1, const std::string &optMethodDeriv=OPTIMIZATION_GRADIENT)
Optimize numerical parameters assuming a global clock (branch heights, substitution model & rate dist...
NNI topology search method.
NaNListener(Optimizer *optimizer, Function *function)