bpp-phyl  2.2.0
bpp::OptimizationTools Class Reference

Optimization methods for phylogenetic inference. More...

#include <Bpp/Phyl/OptimizationTools.h>

+ Collaboration diagram for bpp::OptimizationTools:

Classes

class  ScaleFunction
 

Public Member Functions

 OptimizationTools ()
 
virtual ~OptimizationTools ()
 

Static Public Member Functions

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) throw (Exception)
 Optimize numerical parameters (branch length, substitution model & rate distribution) of a TreeLikelihood function. More...
 
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) throw (Exception)
 Optimize numerical parameters (branch length, substitution model & rate distribution) of a TreeLikelihood function. More...
 
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) throw (Exception)
 Optimize branch lengths parameters of a TreeLikelihood function. More...
 
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) throw (Exception)
 Optimize numerical parameters assuming a global clock (branch heights, substitution model & rate distribution) of a ClockTreeLikelihood function. More...
 
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) throw (Exception)
 Optimize numerical parameters assuming a global clock (branch heights, substitution model & rate distribution) of a ClockTreeLikelihood function. More...
 
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) throw (Exception)
 Optimize the scale of a TreeLikelihood. More...
 
static NNIHomogeneousTreeLikelihoodoptimizeTreeNNI (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) throw (Exception)
 Optimize all parameters from a TreeLikelihood object, including tree topology using Nearest Neighbor Interchanges. More...
 
static NNIHomogeneousTreeLikelihoodoptimizeTreeNNI2 (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) throw (Exception)
 Optimize all parameters from a TreeLikelihood object, including tree topology using Nearest Neighbor Interchanges. More...
 
static DRTreeParsimonyScoreoptimizeTreeNNI (DRTreeParsimonyScore *tp, unsigned int verbose=1)
 Optimize tree topology from a DRTreeParsimonyScore using Nearest Neighbor Interchanges. More...
 
static DistanceMatrix * estimateDistanceMatrix (DistanceEstimation &estimationMethod, const ParameterList &parametersToIgnore, const std::string &param=DISTANCEMETHOD_INIT, unsigned int verbose=0) throw (Exception)
 Estimate a distance matrix using maximum likelihood. More...
 
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) throw (Exception)
 Build a tree using a distance method. More...
 

Static Public Attributes

static std::string OPTIMIZATION_GRADIENT = "gradient"
 
static std::string OPTIMIZATION_NEWTON = "newton"
 
static std::string OPTIMIZATION_BRENT = "Brent"
 
static std::string OPTIMIZATION_BFGS = "BFGS"
 
static std::string DISTANCEMETHOD_INIT = "init"
 
static std::string DISTANCEMETHOD_PAIRWISE = "pairwise"
 
static std::string DISTANCEMETHOD_ITERATIONS = "iterations"
 

Detailed Description

Optimization methods for phylogenetic inference.

This class provides optimization methods for phylogenetics. Parameters of the optimization methods are set to work with TreeLikelihood object. Some non trivial parameters are left to the user choice (tolerance, maximum number of function evaluation, output streams).

Definition at line 300 of file OptimizationTools.h.

Constructor & Destructor Documentation

◆ OptimizationTools()

OptimizationTools::OptimizationTools ( )

Definition at line 64 of file OptimizationTools.cpp.

◆ ~OptimizationTools()

OptimizationTools::~OptimizationTools ( )
virtual

Definition at line 66 of file OptimizationTools.cpp.

Member Function Documentation

◆ buildDistanceTree()

TreeTemplate< Node > * OptimizationTools::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 
)
throw (Exception
)
static

Build a tree using a distance method.

This method estimate a distance matrix using a DistanceEstimation object, and then builds the phylogenetic tree using a AgglomerativeDistanceMethod object. The main issue here is to estimate non-branch lengths parameters, as substitution model and rate distribution parameters. Three options are provideed here:

  • DISTANCEMETHOD_INIT (default) keep parameters to there initial value,
  • DISTANCEMETHOD_PAIRWISE estimated parameters in a pairwise manner, which is standard but not that satisfying...
  • DISTANCEMETHOD_ITERATIONS uses Ninio et al's iterative algorithm, which uses Maximum Likelihood to estimate these parameters, and then update the distance matrix. Ninio M, Privman E, Pupko T, Friedman N. Phylogeny reconstruction: increasing the accuracy of pairwise distance estimation using Bayesian inference of evolutionary rates. Bioinformatics. 2007 Jan 15;23(2):e136-41.
Parameters
estimationMethodThe distance estimation object to use.
reconstructionMethodThe tree reconstruction object to use.
parametersToIgnoreA list of parameters to ignore while optimizing parameters.
optimizeBrLenTell if branch lengths should be optimized together with other parameters. This may lead to more accurate parameter estimation, but is slower.
paramString describing the type of optimization to use.
toleranceThreshold on likelihood for stopping the iterative procedure. Used only with param=DISTANCEMETHOD_ITERATIONS.
tlEvalMaxMaximum number of likelihood computations to perform when optimizing parameters. Used only with param=DISTANCEMETHOD_ITERATIONS.
profilerOutput stream used by optimizer. Used only with param=DISTANCEMETHOD_ITERATIONS.
messengerOutput stream used by optimizer. Used only with param=DISTANCEMETHOD_ITERATIONS.
verboseVerbose level.

Definition at line 677 of file OptimizationTools.cpp.

References DISTANCEMETHOD_ITERATIONS, DISTANCEMETHOD_PAIRWISE, bpp::AbstractHomogeneousTreeLikelihood::initialize(), optimizeNumericalParameters(), and bpp::TreeTools::robinsonFouldsDistance().

◆ estimateDistanceMatrix()

DistanceMatrix * OptimizationTools::estimateDistanceMatrix ( DistanceEstimation estimationMethod,
const ParameterList &  parametersToIgnore,
const std::string &  param = DISTANCEMETHOD_INIT,
unsigned int  verbose = 0 
)
throw (Exception
)
static

Estimate a distance matrix using maximum likelihood.

This method estimate a distance matrix using a DistanceEstimation object. The main issue here is to estimate non-branch lengths parameters, as substitution model and rate distribution parameters. Twoe options are provideed here:

  • DISTANCEMETHOD_INIT (default) keep parameters to there initial value,
  • DISTANCEMETHOD_PAIRWISE estimated parameters in a pairwise manner, which is standard but not that satisfying...
Parameters
estimationMethodThe distance estimation object to use.
parametersToIgnoreA list of parameters to ignore while optimizing parameters.
paramString describing the type of optimization to use.
verboseVerbose level.
See also
buildDistanceTree for a procedure to jointly estimate the distance matrix and underlying tree.

Definition at line 647 of file OptimizationTools.cpp.

References DISTANCEMETHOD_INIT, and DISTANCEMETHOD_PAIRWISE.

◆ optimizeBranchLengthsParameters()

unsigned int OptimizationTools::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 
)
throw (Exception
)
static

Optimize branch lengths parameters of a TreeLikelihood function.

Uses Newton's method.

A condition over function values is used as a stop condition for the algorithm.

See also
NewtonBrentMetaOptimizer
Parameters
tlA pointer toward the TreeLikelihood object to optimize.
parametersThe list of parameters to optimize. The intersection of branch length parameters and the input set will be used. Use tl->getBranchLengthsParameters() in order to estimate all branch length parameters.
listenerA pointer toward an optimization listener, if needed.
toleranceThe tolerance to use in the algorithm.
tlEvalMaxThe maximum number of function evaluations.
messageHandlerThe massage handler.
profilerThe profiler.
verboseThe verbose level.
optMethodDerivOptimization type for derivable parameters (first or second order derivatives).
See also
OPTIMIZATION_NEWTON, OPTIMIZATION_GRADIENT
Exceptions
Exceptionany exception thrown by the Optimizer.

Definition at line 348 of file OptimizationTools.cpp.

References OPTIMIZATION_BFGS, OPTIMIZATION_GRADIENT, and OPTIMIZATION_NEWTON.

◆ optimizeNumericalParameters()

unsigned int OptimizationTools::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 
)
throw (Exception
)
static

Optimize numerical parameters (branch length, substitution model & rate distribution) of a TreeLikelihood function.

Uses Newton's method for branch length and Brent or BFGS one dimensional method for other parameters.

A condition over function values is used as a stop condition for the algorithm.

See also
BrentOneDimension, BFGSMultiDimensions
Parameters
tlA pointer toward the TreeLikelihood object to optimize.
parametersThe list of parameters to optimize. Use tl->getIndependentParameters() in order to estimate all parameters.
listenerA pointer toward an optimization listener, if needed.
nstepThe number of progressive steps to perform (see NewtonBrentMetaOptimizer). 1 means full precision from start.
toleranceThe tolerance to use in the algorithm.
tlEvalMaxThe maximum number of function evaluations.
messageHandlerThe massage handler.
profilerThe profiler.
reparametrizationTell if parameters should be transformed in order to remove constraints. This can improve optimization, but is a bit slower.
verboseThe verbose level.
optMethodDerivOptimization type for derivable parameters (first or second order derivatives).
See also
OPTIMIZATION_NEWTON, OPTIMIZATION_GRADIENT
Parameters
optMethodModelOptimization type for model parameters (Brent or BFGS).
See also
OPTIMIZATION_BRENT, OPTIMIZATION_BFGS
Exceptions
Exceptionany exception thrown by the Optimizer.

Definition at line 150 of file OptimizationTools.cpp.

References OPTIMIZATION_BFGS, OPTIMIZATION_BRENT, OPTIMIZATION_GRADIENT, and OPTIMIZATION_NEWTON.

Referenced by buildDistanceTree(), bpp::PhylogeneticsApplicationTools::optimizeParameters(), optimizeTreeNNI(), and bpp::NNITopologyListener::topologyChangeSuccessful().

◆ optimizeNumericalParameters2()

unsigned int OptimizationTools::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 
)
throw (Exception
)
static

Optimize numerical parameters (branch length, substitution model & rate distribution) of a TreeLikelihood function.

Uses Newton's method for all parameters, branch length derivatives are computed analytically, derivatives for other parameters numerically.

See also
PseudoNewtonOptimizer
Parameters
tlA pointer toward the TreeLikelihood object to optimize.
parametersThe list of parameters to optimize. Use tl->getIndependentParameters() in order to estimate all parameters.
listenerA pointer toward an optimization listener, if needed.
toleranceThe tolerance to use in the algorithm.
tlEvalMaxThe maximum number of function evaluations.
messageHandlerThe massage handler.
profilerThe profiler.
reparametrizationTell if parameters should be transformed in order to remove constraints. This can improve optimization, but is a bit slower.
useClockTell if branch lengths have to be optimized under a global molecular clock constraint.
verboseThe verbose level.
optMethodDerivOptimization type for derivable parameters (first or second order derivatives).
See also
OPTIMIZATION_NEWTON, OPTIMIZATION_GRADIENT
Exceptions
Exceptionany exception thrown by the Optimizer.

Definition at line 255 of file OptimizationTools.cpp.

References bpp::GlobalClockTreeLikelihoodFunctionWrapper::getHeightParameters(), OPTIMIZATION_BFGS, OPTIMIZATION_GRADIENT, and OPTIMIZATION_NEWTON.

Referenced by bpp::PhylogeneticsApplicationTools::optimizeParameters(), optimizeTreeNNI2(), and bpp::NNITopologyListener2::topologyChangeSuccessful().

◆ optimizeNumericalParametersWithGlobalClock()

unsigned int OptimizationTools::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 
)
throw (Exception
)
static

Optimize numerical parameters assuming a global clock (branch heights, substitution model & rate distribution) of a ClockTreeLikelihood function.

Uses Newton or conjugate gradient method for branch length and Brent's one dimensional method for other parameters (NewtonBrentMetaOptimizer). Derivatives are computed analytically.

A condition over function values is used as a stop condition for the algorithm.

See also
NewtonBrentMetaOptimizer
Parameters
clA pointer toward the ClockTreeLikelihood object to optimize.
parametersThe list of parameters to optimize. Use cl->getIndependentParameters() in order to estimate all parameters.
listenerA pointer toward an optimization listener, if needed.
nstepThe number of progressive steps to perform (see NewtonBrentMetaOptimizer). 1 means full precision from start.
toleranceThe tolerance to use in the algorithm.
tlEvalMaxThe maximum number of function evaluations.
messageHandlerThe massage handler.
profilerThe profiler.
verboseThe verbose level.
optMethodDerivOptimization type for derivable parameters (first or second order derivatives).
See also
OPTIMIZATION_NEWTON, OPTIMIZATION_GRADIENT
Exceptions
Exceptionany exception thrown by the Optimizer.
Deprecated:
See optimizeNumericalParameters2 as a more general replacement.

Definition at line 406 of file OptimizationTools.cpp.

References OPTIMIZATION_GRADIENT, and OPTIMIZATION_NEWTON.

Referenced by bpp::PhylogeneticsApplicationTools::optimizeParameters().

◆ optimizeNumericalParametersWithGlobalClock2()

unsigned int OptimizationTools::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 
)
throw (Exception
)
static

Optimize numerical parameters assuming a global clock (branch heights, substitution model & rate distribution) of a ClockTreeLikelihood function.

Uses Newton or conjugate gradient method for all parameters, branch length derivatives are computed analytically, derivatives for other parameters numerically.

See also
PseudoNewtonOptimizer
Parameters
clA pointer toward the ClockTreeLikelihood object to optimize.
parametersThe list of parameters to optimize. Use cl->getIndependentParameters() in order to estimate all parameters.
listenerA pointer toward an optimization listener, if needed.
toleranceThe tolerance to use in the algorithm.
tlEvalMaxThe maximum number of function evaluations.
messageHandlerThe massage handler.
profilerThe profiler.
verboseThe verbose level.
optMethodDerivOptimization type for derivable parameters (first or second order derivatives).
See also
OPTIMIZATION_NEWTON, OPTIMIZATION_GRADIENT
Exceptions
Exceptionany exception thrown by the Optimizer.
Deprecated:
See optimizeNumericalParameters2 as a more general replacement.

Definition at line 476 of file OptimizationTools.cpp.

References OPTIMIZATION_GRADIENT, and OPTIMIZATION_NEWTON.

Referenced by bpp::PhylogeneticsApplicationTools::optimizeParameters().

◆ optimizeTreeNNI() [1/2]

NNIHomogeneousTreeLikelihood * OptimizationTools::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 
)
throw (Exception
)
static

Optimize all parameters from a TreeLikelihood object, including tree topology using Nearest Neighbor Interchanges.

This function takes as input a TreeLikelihood object implementing the NNISearchable interface.

Details: A NNITopologySearch object is instanciated and is associated an additional TopologyListener. This listener is used to re-estimate numerical parameters after one or several topology change. By default, the PHYML option is used for the NNITopologySearch object, and numerical parameters are re-estimated every 4 NNI runs (as in the phyml software).

The optimizeNumericalParameters method is used for estimating numerical parameters. The tolerance passed to this function is specified as input parameters. They are generally very high to avoid local optima.

Parameters
tlA pointer toward the TreeLikelihood object to optimize.
parametersThe list of parameters to optimize. Use tl->getIndependentParameters() in order to estimate all parameters.
optimizeNumFirstTell if we must optimize numerical parameters before searching topology.
tolBeforeThe tolerance to use when estimating numerical parameters before topology search (if optimizeNumFirst is set to 'true').
tolDuringThe tolerance to use when estimating numerical parameters during the topology search.
tlEvalMaxThe maximum number of function evaluations.
numStepNumber of NNI rounds before re-estimating numerical parameters.
messageHandlerThe massage handler.
profilerThe profiler.
reparametrizationTell if parameters should be transformed in order to remove constraints. This can improve optimization, but is a bit slower.
verboseThe verbose level.
optMethodOption passed to optimizeNumericalParameters.
nStepOption passed to optimizeNumericalParameters.
nniMethodNNI algorithm to use.
Returns
A pointer toward the final likelihood object. This pointer may be the same as passed in argument (tl), but in some cases the algorithm clone this object. We may change this bahavior in the future... You hence should write something like
Exceptions
Exceptionany exception thrown by the optimizer.

Definition at line 564 of file OptimizationTools.cpp.

References bpp::NNITopologySearch::addTopologyListener(), bpp::NNITopologySearch::getSearchableObject(), optimizeNumericalParameters(), bpp::NNITopologySearch::search(), and bpp::NNITopologyListener::setNumericalOptimizationCounter().

Referenced by bpp::TreeTools::MRP(), bpp::TreeTools::MRPMultilabel(), and bpp::PhylogeneticsApplicationTools::optimizeParameters().

◆ optimizeTreeNNI() [2/2]

DRTreeParsimonyScore * OptimizationTools::optimizeTreeNNI ( DRTreeParsimonyScore tp,
unsigned int  verbose = 1 
)
static

Optimize tree topology from a DRTreeParsimonyScore using Nearest Neighbor Interchanges.

Parameters
tpA pointer toward the DRTreeParsimonyScore object to optimize.
verboseThe verbose level.
Returns
A pointer toward the final parsimony score object. This pointer may be the same as passed in argument (tl), but in some cases the algorithm clone this object. We may change this bahavior in the future... You hence should write something like

Definition at line 629 of file OptimizationTools.cpp.

References bpp::NNITopologySearch::getSearchableObject(), bpp::NNITopologySearch::PHYML, and bpp::NNITopologySearch::search().

◆ optimizeTreeNNI2()

NNIHomogeneousTreeLikelihood * OptimizationTools::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 
)
throw (Exception
)
static

Optimize all parameters from a TreeLikelihood object, including tree topology using Nearest Neighbor Interchanges.

This function takes as input a TreeLikelihood object implementing the NNISearchable interface.

Details: A NNITopologySearch object is instanciated and is associated an additional TopologyListener. This listener is used to re-estimate numerical parameters after one or several topology change. By default, the PHYML option is used for the NNITopologySearch object, and numerical parameters are re-estimated every 4 NNI runs (as in the phyml software).

The optimizeNumericalParameters2 method is used for estimating numerical parameters. The tolerance passed to this function is specified as input parameters. They are generally very high to avoid local optima.

Parameters
tlA pointer toward the TreeLikelihood object to optimize.
parametersThe list of parameters to optimize. Use tl->getIndependentParameters() in order to estimate all parameters.
optimizeNumFirstTell if we must optimize numerical parameters before searching topology.
tolBeforeThe tolerance to use when estimating numerical parameters before topology search (if optimizeNumFirst is set to 'true').
tolDuringThe tolerance to use when estimating numerical parameters during the topology search.
tlEvalMaxThe maximum number of function evaluations.
numStepNumber of NNI rounds before re-estimating numerical parameters.
messageHandlerThe massage handler.
profilerThe profiler.
reparametrizationTell if parameters should be transformed in order to remove constraints. This can improve optimization, but is a bit slower.
verboseThe verbose level.
optMethodOption passed to optimizeNumericalParameters2.
nniMethodNNI algorithm to use.
Returns
A pointer toward the final likelihood object. This pointer may be the same as passed in argument (tl), but in some cases the algorithm clone this object. We may change this bahavior in the future... You hence should write something like
Exceptions
Exceptionany exception thrown by the optimizer.

Definition at line 597 of file OptimizationTools.cpp.

References bpp::NNITopologySearch::addTopologyListener(), bpp::NNITopologySearch::getSearchableObject(), optimizeNumericalParameters2(), bpp::NNITopologySearch::search(), and bpp::NNITopologyListener2::setNumericalOptimizationCounter().

Referenced by bpp::PhylogeneticsApplicationTools::optimizeParameters().

◆ optimizeTreeScale()

unsigned int OptimizationTools::optimizeTreeScale ( TreeLikelihood tl,
double  tolerance = 0.000001,
unsigned int  tlEvalMax = 1000000,
OutputStream *  messageHandler = ApplicationTools::message,
OutputStream *  profiler = ApplicationTools::message,
unsigned int  verbose = 1 
)
throw (Exception
)
static

Optimize the scale of a TreeLikelihood.

This method only works on branch lengths parameters. It multiply all branch length by a factor 'x' which is optimized using Brent's algorithm in one dimension. This method may be usefull for scaling a tree whose branch lengths come from the Neighbor-Joining algorithm for instance.

Practically, and contrarily to what one may expect, this does not speed up the optimization!

A condition over parameters is used as a stop condition for the algorithm.

Parameters
tlA pointer toward the TreeLikelihood object to optimize.
toleranceThe tolerance to use in the algorithm.
tlEvalMaxThe maximum number of function evaluations.
messageHandlerThe massage handler.
profilerThe profiler.
verboseThe verbose level.
Exceptions
Exceptionany exception thrown by the optimizer.

Definition at line 121 of file OptimizationTools.cpp.

References bpp::OptimizationTools::ScaleFunction::getParameters().

Referenced by bpp::PhylogeneticsApplicationTools::optimizeParameters().

Member Data Documentation

◆ DISTANCEMETHOD_INIT

std::string OptimizationTools::DISTANCEMETHOD_INIT = "init"
static

Definition at line 763 of file OptimizationTools.h.

Referenced by estimateDistanceMatrix().

◆ DISTANCEMETHOD_ITERATIONS

std::string OptimizationTools::DISTANCEMETHOD_ITERATIONS = "iterations"
static

Definition at line 765 of file OptimizationTools.h.

Referenced by buildDistanceTree().

◆ DISTANCEMETHOD_PAIRWISE

std::string OptimizationTools::DISTANCEMETHOD_PAIRWISE = "pairwise"
static

Definition at line 764 of file OptimizationTools.h.

Referenced by buildDistanceTree(), and estimateDistanceMatrix().

◆ OPTIMIZATION_BFGS

std::string OptimizationTools::OPTIMIZATION_BFGS = "BFGS"
static

◆ OPTIMIZATION_BRENT

std::string OptimizationTools::OPTIMIZATION_BRENT = "Brent"
static

◆ OPTIMIZATION_GRADIENT

◆ OPTIMIZATION_NEWTON


The documentation for this class was generated from the following files: