bpp-phyl3  3.0.0
OptimizationTools.h
Go to the documentation of this file.
1 // SPDX-FileCopyrightText: The Bio++ Development Group
2 //
3 // SPDX-License-Identifier: CECILL-2.1
4 
5 #ifndef BPP_PHYL_LEGACY_OPTIMIZATIONTOOLS_H
6 #define BPP_PHYL_LEGACY_OPTIMIZATIONTOOLS_H
7 
9 #include <Bpp/Io/OutputStream.h>
12 
13 #include "../Parsimony/DRTreeParsimonyScore.h"
14 #include "../Tree/TreeTemplate.h"
16 #include "Tree/NNITopologySearch.h"
17 
18 namespace bpp
19 {
24  public virtual TopologyListener
25 {
26 private:
27  std::shared_ptr<NNITopologySearch> topoSearch_;
29  double tolerance_;
30  std::shared_ptr<OutputStream> messenger_;
31  std::shared_ptr<OutputStream> profiler_;
32  unsigned int verbose_;
33  unsigned int optimizeCounter_;
34  unsigned int optimizeNumerical_;
35  std::string optMethod_;
36  unsigned int nStep_;
38 
39 public:
58  std::shared_ptr<NNITopologySearch> ts,
59  const ParameterList& parameters,
60  double tolerance,
61  std::shared_ptr<OutputStream> messenger,
62  std::shared_ptr<OutputStream> profiler,
63  unsigned int verbose,
64  const std::string& optMethod,
65  unsigned int nStep,
66  bool reparametrization) :
67  topoSearch_(ts),
68  parameters_(parameters),
69  tolerance_(tolerance),
70  messenger_(messenger),
71  profiler_(profiler),
72  verbose_(verbose),
75  optMethod_(optMethod),
76  nStep_(nStep),
77  reparametrization_(reparametrization) {}
78 
84  profiler_(tl.profiler_),
85  verbose_(tl.verbose_),
89  nStep_(tl.nStep_),
91  {}
92 
94  {
99  profiler_ = tl.profiler_;
100  verbose_ = tl.verbose_;
103  optMethod_ = tl.optMethod_;
104  nStep_ = tl.nStep_;
106  return *this;
107  }
108 
109  NNITopologyListener* clone() const { return new NNITopologyListener(*this); }
110 
111  virtual ~NNITopologyListener() {}
112 
113 public:
117 };
118 
123  public TopologyListener
124 {
125 private:
126  std::shared_ptr<NNITopologySearch> topoSearch_;
128  double tolerance_;
129  std::shared_ptr<OutputStream> messenger_;
130  std::shared_ptr<OutputStream> profiler_;
131  unsigned int verbose_;
132  unsigned int optimizeCounter_;
133  unsigned int optimizeNumerical_;
134  std::string optMethod_;
136 
137 public:
155  std::shared_ptr<NNITopologySearch> ts,
156  const ParameterList& parameters,
157  double tolerance,
158  std::shared_ptr<OutputStream> messenger,
159  std::shared_ptr<OutputStream> profiler,
160  unsigned int verbose,
161  const std::string& optMethod,
162  bool reparametrization) :
163  topoSearch_(ts),
164  parameters_(parameters),
165  tolerance_(tolerance),
166  messenger_(messenger),
167  profiler_(profiler),
168  verbose_(verbose),
169  optimizeCounter_(0),
171  optMethod_(optMethod),
172  reparametrization_(reparametrization) {}
173 
179  profiler_(tl.profiler_),
180  verbose_(tl.verbose_),
185  {}
186 
188  {
191  tolerance_ = tl.tolerance_;
192  messenger_ = tl.messenger_;
193  profiler_ = tl.profiler_;
194  verbose_ = tl.verbose_;
197  optMethod_ = tl.optMethod_;
199  return *this;
200  }
201 
202  NNITopologyListener2* clone() const { return new NNITopologyListener2(*this); }
203 
205 
206 public:
210 };
211 
212 
222  public OptimizationTools
223 {
224 public:
226  virtual ~LegacyOptimizationTools();
227 
228 public:
255  static unsigned int optimizeNumericalParameters(
256  std::shared_ptr<DiscreteRatesAcrossSitesTreeLikelihoodInterface> tl,
257  const ParameterList& parameters,
258  std::shared_ptr<OptimizationListener> listener = nullptr,
259  unsigned int nstep = 1,
260  double tolerance = 0.000001,
261  unsigned int tlEvalMax = 1000000,
262  std::shared_ptr<OutputStream> messageHandler = ApplicationTools::message,
263  std::shared_ptr<OutputStream> profiler = ApplicationTools::message,
264  bool reparametrization = false,
265  unsigned int verbose = 1,
266  const std::string& optMethodDeriv = OPTIMIZATION_NEWTON,
267  const std::string& optMethodModel = OPTIMIZATION_BRENT);
268 
291  static unsigned int optimizeNumericalParameters2(
292  std::shared_ptr<DiscreteRatesAcrossSitesTreeLikelihoodInterface> tl,
293  const ParameterList& parameters,
294  std::shared_ptr<OptimizationListener> listener = nullptr,
295  double tolerance = 0.000001,
296  unsigned int tlEvalMax = 1000000,
297  std::shared_ptr<OutputStream> messageHandler = ApplicationTools::message,
298  std::shared_ptr<OutputStream> profiler = ApplicationTools::message,
299  bool reparametrization = false,
300  bool useClock = false,
301  unsigned int verbose = 1,
302  const std::string& optMethodDeriv = OPTIMIZATION_NEWTON);
303 
325  static unsigned int optimizeBranchLengthsParameters(
326  std::shared_ptr<DiscreteRatesAcrossSitesTreeLikelihoodInterface> tl,
327  const ParameterList& parameters,
328  std::shared_ptr<OptimizationListener> listener = nullptr,
329  double tolerance = 0.000001,
330  unsigned int tlEvalMax = 1000000,
331  std::shared_ptr<OutputStream> messageHandler = ApplicationTools::message,
332  std::shared_ptr<OutputStream> profiler = ApplicationTools::message,
333  unsigned int verbose = 1,
334  const std::string& optMethodDeriv = OPTIMIZATION_NEWTON);
335 
336 private:
338  public virtual FunctionInterface,
339  public ParametrizableAdapter
340  {
341 private:
342  std::shared_ptr<TreeLikelihoodInterface> tl_;
344 
345 public:
346  ScaleFunction(std::shared_ptr<TreeLikelihoodInterface> tl);
347 
349  tl_(sf.tl_),
350  brLen_(sf.brLen_),
351  lambda_(sf.lambda_)
352  {}
353 
355  {
356  tl_ = sf.tl_;
357  brLen_ = sf.brLen_;
358  lambda_ = sf.lambda_;
359  return *this;
360  }
361 
362  virtual ~ScaleFunction();
363 
364  ScaleFunction* clone() const { return new ScaleFunction(*this); }
365 
366 public:
367  void setParameters(const ParameterList& lambda);
368  double getValue() const;
369  const ParameterList& getParameters() const { return lambda_; }
370  const Parameter& parameter(const std::string& name) const
371  {
372  if (name == "lambda") return lambda_[0];
373  else throw ParameterNotFoundException("ScaleFunction::getParameter.", name);
374  }
375  double getParameterValue(const std::string& name) const
376  {
377  return lambda_.parameter(name).getValue();
378  }
379  size_t getNumberOfParameters() const { return 1; }
380  size_t getNumberOfIndependentParameters() const { return 1; }
381 
382 protected:
384  };
385 
386 public:
409  static unsigned int optimizeTreeScale(
410  std::shared_ptr<TreeLikelihoodInterface> tl,
411  double tolerance = 0.000001,
412  unsigned int tlEvalMax = 1000000,
413  std::shared_ptr<OutputStream> messageHandler = ApplicationTools::message,
414  std::shared_ptr<OutputStream> profiler = ApplicationTools::message,
415  unsigned int verbose = 1);
416 
417 
457  static std::shared_ptr<NNIHomogeneousTreeLikelihood> optimizeTreeNNI(
458  std::shared_ptr<NNIHomogeneousTreeLikelihood> tl,
459  const ParameterList& parameters,
460  bool optimizeNumFirst = true,
461  double tolBefore = 100,
462  double tolDuring = 100,
463  unsigned int tlEvalMax = 1000000,
464  unsigned int numStep = 1,
465  std::shared_ptr<OutputStream> messageHandler = ApplicationTools::message,
466  std::shared_ptr<OutputStream> profiler = ApplicationTools::message,
467  bool reparametrization = false,
468  unsigned int verbose = 1,
469  const std::string& optMethod = OptimizationTools::OPTIMIZATION_NEWTON,
470  unsigned int nStep = 1,
471  const std::string& nniMethod = NNITopologySearch::PHYML);
472 
511  static std::shared_ptr<NNIHomogeneousTreeLikelihood> optimizeTreeNNI2(
512  std::shared_ptr<NNIHomogeneousTreeLikelihood> tl,
513  const ParameterList& parameters,
514  bool optimizeNumFirst = true,
515  double tolBefore = 100,
516  double tolDuring = 100,
517  unsigned int tlEvalMax = 1000000,
518  unsigned int numStep = 1,
519  std::shared_ptr<OutputStream> messageHandler = ApplicationTools::message,
520  std::shared_ptr<OutputStream> profiler = ApplicationTools::message,
521  bool reparametrization = false,
522  unsigned int verbose = 1,
523  const std::string& optMethod = OptimizationTools::OPTIMIZATION_NEWTON,
524  const std::string& nniMethod = NNITopologySearch::PHYML);
525 
539  static std::shared_ptr<DRTreeParsimonyScore> optimizeTreeNNI(
540  std::shared_ptr<DRTreeParsimonyScore> tp,
541  unsigned int verbose = 1);
542 };
543 } // end of namespace bpp.
544 #endif // BPP_PHYL_LEGACY_OPTIMIZATIONTOOLS_H
static std::shared_ptr< OutputStream > message
void setParameters(const ParameterList &lambda)
const Parameter & parameter(const std::string &name) const
ScaleFunction(std::shared_ptr< TreeLikelihoodInterface > tl)
std::shared_ptr< TreeLikelihoodInterface > tl_
const ParameterList & getParameters() const
double getParameterValue(const std::string &name) const
ScaleFunction & operator=(const ScaleFunction &sf)
Optimization methods for phylogenetic inference.
static std::shared_ptr< NNIHomogeneousTreeLikelihood > optimizeTreeNNI2(std::shared_ptr< NNIHomogeneousTreeLikelihood > tl, const ParameterList &parameters, bool optimizeNumFirst=true, double tolBefore=100, double tolDuring=100, unsigned int tlEvalMax=1000000, unsigned int numStep=1, std::shared_ptr< OutputStream > messageHandler=ApplicationTools::message, std::shared_ptr< 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 ...
static unsigned int optimizeNumericalParameters2(std::shared_ptr< DiscreteRatesAcrossSitesTreeLikelihoodInterface > tl, const ParameterList &parameters, std::shared_ptr< OptimizationListener > listener=nullptr, double tolerance=0.000001, unsigned int tlEvalMax=1000000, std::shared_ptr< OutputStream > messageHandler=ApplicationTools::message, std::shared_ptr< 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...
static unsigned int optimizeTreeScale(std::shared_ptr< TreeLikelihoodInterface > tl, double tolerance=0.000001, unsigned int tlEvalMax=1000000, std::shared_ptr< OutputStream > messageHandler=ApplicationTools::message, std::shared_ptr< OutputStream > profiler=ApplicationTools::message, unsigned int verbose=1)
Optimize the scale of a TreeLikelihood.
static std::shared_ptr< NNIHomogeneousTreeLikelihood > optimizeTreeNNI(std::shared_ptr< NNIHomogeneousTreeLikelihood > tl, const ParameterList &parameters, bool optimizeNumFirst=true, double tolBefore=100, double tolDuring=100, unsigned int tlEvalMax=1000000, unsigned int numStep=1, std::shared_ptr< OutputStream > messageHandler=ApplicationTools::message, std::shared_ptr< 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 ...
static unsigned int optimizeNumericalParameters(std::shared_ptr< DiscreteRatesAcrossSitesTreeLikelihoodInterface > tl, const ParameterList &parameters, std::shared_ptr< OptimizationListener > listener=nullptr, unsigned int nstep=1, double tolerance=0.000001, unsigned int tlEvalMax=1000000, std::shared_ptr< OutputStream > messageHandler=ApplicationTools::message, std::shared_ptr< 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...
static unsigned int optimizeBranchLengthsParameters(std::shared_ptr< DiscreteRatesAcrossSitesTreeLikelihoodInterface > tl, const ParameterList &parameters, std::shared_ptr< OptimizationListener > listener=nullptr, double tolerance=0.000001, unsigned int tlEvalMax=1000000, std::shared_ptr< OutputStream > messageHandler=ApplicationTools::message, std::shared_ptr< OutputStream > profiler=ApplicationTools::message, unsigned int verbose=1, const std::string &optMethodDeriv=OPTIMIZATION_NEWTON)
Optimize branch lengths parameters of a TreeLikelihood function.
Listener used internally by the optimizeTreeNNI2 method.
NNITopologyListener2 & operator=(const NNITopologyListener2 &tl)
NNITopologyListener2(const NNITopologyListener2 &tl)
std::shared_ptr< OutputStream > messenger_
NNITopologyListener2(std::shared_ptr< NNITopologySearch > ts, const ParameterList &parameters, double tolerance, std::shared_ptr< OutputStream > messenger, std::shared_ptr< OutputStream > profiler, unsigned int verbose, const std::string &optMethod, bool reparametrization)
Build a new NNITopologyListener2 object.
void topologyChangeTested(const TopologyChangeEvent &event)
Notify a topology change event.
std::shared_ptr< OutputStream > profiler_
std::shared_ptr< NNITopologySearch > topoSearch_
void setNumericalOptimizationCounter(unsigned int c)
void topologyChangeSuccessful(const TopologyChangeEvent &event)
Tell that a topology change is definitive.
NNITopologyListener2 * clone() const
Listener used internally by the optimizeTreeNNI method.
std::shared_ptr< OutputStream > messenger_
NNITopologyListener * clone() const
void setNumericalOptimizationCounter(unsigned int c)
NNITopologyListener(std::shared_ptr< NNITopologySearch > ts, const ParameterList &parameters, double tolerance, std::shared_ptr< OutputStream > messenger, std::shared_ptr< OutputStream > profiler, unsigned int verbose, const std::string &optMethod, unsigned int nStep, bool reparametrization)
Build a new NNITopologyListener object.
std::shared_ptr< OutputStream > profiler_
void topologyChangeSuccessful(const TopologyChangeEvent &event)
Tell that a topology change is definitive.
NNITopologyListener(const NNITopologyListener &tl)
NNITopologyListener & operator=(const NNITopologyListener &tl)
std::shared_ptr< NNITopologySearch > topoSearch_
void topologyChangeTested(const TopologyChangeEvent &event)
Notify a topology change event.
static const std::string PHYML
Optimization methods for phylogenetic inference.
static std::string OPTIMIZATION_NEWTON
static std::string OPTIMIZATION_BRENT
virtual const Parameter & parameter(const std::string &name) const
virtual double getValue() const
Class for notifying new toplogy change events.
Implement this interface to be notified when the topology of a tree has changed during topology searc...
Defines the basic types of data flow nodes.