bpp-phyl  2.2.0
bpp::TwoTreeLikelihood Class Reference

This class is a simplified version of DRHomogeneousTreeLikelihood for 2-Trees. More...

#include <Bpp/Phyl/Distance/DistanceEstimation.h>

+ Inheritance diagram for bpp::TwoTreeLikelihood:
+ Collaboration diagram for bpp::TwoTreeLikelihood:

Public Member Functions

 TwoTreeLikelihood (const std::string &seq1, const std::string &seq2, const SiteContainer &data, SubstitutionModel *model, DiscreteDistribution *rDist, bool verbose) throw (Exception)
 
 TwoTreeLikelihood (const TwoTreeLikelihood &lik)
 
TwoTreeLikelihoodoperator= (const TwoTreeLikelihood &lik)
 
TwoTreeLikelihoodclone () const
 
virtual ~TwoTreeLikelihood ()
 
const SubstitutionModelgetSubstitutionModel () const
 Get the substitution model used for the computation. More...
 
SubstitutionModelgetSubstitutionModel ()
 Get the substitution model used for the computation. More...
 
ConstBranchModelIteratorgetNewBranchModelIterator (int nodeId) const throw (NotImplementedException)
 
ConstSiteModelIteratorgetNewSiteModelIterator (size_t siteIndex) const throw (NotImplementedException)
 
void setParameters (const ParameterList &parameters) throw (ParameterNotFoundException, ConstraintException)
 Implements the Function interface. More...
 
double getValue () const throw (Exception)
 
virtual void initBranchLengthsParameters ()
 
virtual void setMinimumBranchLength (double minimum)
 
virtual double getMinimumBranchLength () const
 
The TreeLikelihood interface.

Other methods are implemented in the AbstractTreeLikelihood class.

size_t getNumberOfStates () const
 
const std::vector< int > & getAlphabetStates () const
 
int getAlphabetStateAsInt (size_t i) const
 
std::string getAlphabetStateAsChar (size_t i) const
 
TreeLikelihoodDatagetLikelihoodData () throw (NotImplementedException)
 
const TreeLikelihoodDatagetLikelihoodData () const throw (NotImplementedException)
 
double getLikelihood () const
 Get the likelihood for the whole dataset. More...
 
double getLogLikelihood () const
 Get the logarithm of the likelihood for the whole dataset. More...
 
double getLikelihoodForASite (size_t site) const
 Get the likelihood for a site. More...
 
double getLogLikelihoodForASite (size_t site) const
 Get the logarithm of the likelihood for a site. More...
 
ParameterList getBranchLengthsParameters () const
 Get the branch lengths parameters. More...
 
ParameterList getSubstitutionModelParameters () const
 Get the parameters associated to substitution model(s). More...
 
SubstitutionModelgetSubstitutionModel (int nodeId, size_t siteIndex) throw (NodeNotFoundException)
 Get the substitution model associated to a given node and alignment column. More...
 
const SubstitutionModelgetSubstitutionModel (int nodeId, size_t siteIndex) const throw (NodeNotFoundException)
 Get the substitution model associated to a given node and alignment column. More...
 
const std::vector< double > & getRootFrequencies (size_t siteIndex) const
 Get the values of the frequencies for each state in the alphabet at the root node. More...
 
size_t getSiteIndex (size_t site) const throw (IndexOutOfBoundsException)
 Get the index (used for inner computations) of a given site (original alignment column). More...
 
VVVdouble getTransitionProbabilitiesPerRateClass (int nodeId, size_t siteIndex) const
 This method is not applicable for this object. More...
 
void setData (const SiteContainer &sites) throw (Exception)
 Set the dataset for which the likelihood must be evaluated. More...
 
void initialize () throw (Exception)
 Init the likelihood object. More...
 
The DiscreteRatesAcrossSites interface implementation:
double getLikelihoodForASiteForARateClass (size_t site, size_t rateClass) const
 Get the likelihood for a site knowing its rate class. More...
 
double getLogLikelihoodForASiteForARateClass (size_t site, size_t rateClass) const
 Get the logarithm of the likelihood for a site knowing its rate class. More...
 
double getLikelihoodForASiteForARateClassForAState (size_t site, size_t rateClass, int state) const
 Get the likelihood for a site knowing its rate class and its ancestral state. More...
 
double getLogLikelihoodForASiteForARateClassForAState (size_t site, size_t rateClass, int state) const
 Get the logarithm of the likelihood for a site knowing its rate class and its ancestral state. More...
 
DerivableFirstOrder interface.
double getFirstOrderDerivative (const std::string &variable) const throw (Exception)
 
DerivableSecondOrder interface.
double getSecondOrderDerivative (const std::string &variable) const throw (Exception)
 
double getSecondOrderDerivative (const std::string &variable1, const std::string &variable2) const throw (Exception)
 
The TreeLikelihood interface.

Other methods are implemented in the AbstractTreeLikelihood class.

double getLikelihoodForASiteForAState (size_t site, int state) const
 Get the likelihood for a site and for a state. More...
 
double getLogLikelihoodForASiteForAState (size_t site, int state) const
 Get the logarithm of the likelihood for a site and for a state. More...
 
ParameterList getDerivableParameters () const
 All derivable parameters. More...
 
ParameterList getNonDerivableParameters () const
 All non derivable parameters. More...
 
VVdouble getTransitionProbabilities (int nodeId, size_t siteIndex) const
 Retrieves all Pij(t) for a particular branch, defined by the upper node and site. More...
 
The DiscreteRatesAcrossSites interface implementation:
const DiscreteDistribution * getRateDistribution () const
 Get the rate distribution used for the computation. More...
 
DiscreteDistribution * getRateDistribution ()
 Get the rate distribution used for the computation. More...
 
size_t getNumberOfClasses () const
 Get the number of classes. More...
 
ParameterList getRateDistributionParameters () const
 Get the parameters associated to the rate distirbution. More...
 
VVdouble getLikelihoodForEachSiteForEachRateClass () const
 Get the likelihood for each site and each rate class. More...
 
VVdouble getLogLikelihoodForEachSiteForEachRateClass () const
 Get the logarithm of the likelihood for each site and each rate class. More...
 
VVVdouble getLikelihoodForEachSiteForEachRateClassForEachState () const
 Get the likelihood for each site and each rate class and each state. More...
 
VVVdouble getLogLikelihoodForEachSiteForEachRateClassForEachState () const
 Get the logarithm of the likelihood for each site and each rate class and each state. More...
 
VVdouble getPosteriorProbabilitiesOfEachRate () const
 Get the posterior probability for each site of belonging to a particular rate class. More...
 
Vdouble getRateWithMaxPostProbOfEachSite () const
 Get the posterior rate (the one with maximum posterior probability) for each site. More...
 
std::vector< size_t > getRateClassWithMaxPostProbOfEachSite () const
 Get the posterior rate class (the one with maximum posterior probability) for each site. More...
 
Vdouble getPosteriorRateOfEachSite () const
 Get the posterior rate, i.e. averaged over all classes and weighted with posterior probabilities, for each site. More...
 
The TreeLikelihood interface.
const SiteContainer * getData () const
 Get the dataset for which the likelihood must be evaluated. More...
 
const Alphabet * getAlphabet () const
 Get the alphabet associated to the dataset. More...
 
Vdouble getLikelihoodForEachSite () const
 Get the likelihood for each site. More...
 
Vdouble getLogLikelihoodForEachSite () const
 Get the logarithm of the likelihood for each site. More...
 
VVdouble getLikelihoodForEachSiteForEachState () const
 Get the likelihood for each site and for each state. More...
 
VVdouble getLogLikelihoodForEachSiteForEachState () const
 Get the logarithm of the likelihood for each site and for each state. More...
 
size_t getNumberOfSites () const
 Get the number of sites in the dataset. More...
 
const TreegetTree () const
 Get the tree (topology and branch lengths). More...
 
void enableDerivatives (bool yn)
 Tell if derivatives must be computed. More...
 
void enableFirstOrderDerivatives (bool yn)
 
bool enableFirstOrderDerivatives () const
 
void enableSecondOrderDerivatives (bool yn)
 
bool enableSecondOrderDerivatives () const
 
bool isInitialized () const
 

Static Public Member Functions

Generic tools to deal with likelihood arrays
static void resetLikelihoodArray (VVVdouble &likelihoodArray)
 Set all conditional likelihoods to 1. More...
 
static void displayLikelihoodArray (const VVVdouble &likelihoodArray)
 Print the likelihood array to terminal (debugging tool). More...
 

Protected Member Functions

virtual void initTreeLikelihoods (const SequenceContainer &sequences) throw (Exception)
 This method initializes the leaves according to a sequence container. More...
 
void fireParameterChanged (const ParameterList &params)
 
virtual void computeTreeLikelihood ()
 
virtual void computeTreeDLikelihood ()
 
virtual void computeTreeD2Likelihood ()
 
virtual void initParameters ()
 This builds the parameters list from all parametrizable objects, i.e. substitution model, rate distribution and tree. More...
 
virtual void applyParameters () throw (Exception)
 All parameters are stores in a parameter list. More...
 

Protected Attributes

DiscreteDistribution * rateDistribution_
 
const SiteContainer * data_
 
TreeTemplate< Node > * tree_
 
bool computeFirstOrderDerivatives_
 
bool computeSecondOrderDerivatives_
 
bool initialized_
 

Private Attributes

SiteContainer * shrunkData_
 
std::vector< std::string > seqnames_
 
SubstitutionModelmodel_
 
ParameterList brLenParameters_
 
VVVdouble pxy_
 
VVVdouble dpxy_
 
VVVdouble d2pxy_
 
std::vector< size_t > rootPatternLinks_
 As previous, but for the global container. More...
 
std::vector< unsigned int > rootWeights_
 The frequency of each site. More...
 
size_t nbSites_
 
size_t nbClasses_
 
size_t nbStates_
 
size_t nbDistinctSites_
 
VVVdouble rootLikelihoods_
 
VVdouble rootLikelihoodsS_
 
Vdouble rootLikelihoodsSR_
 
Vdouble dLikelihoods_
 
Vdouble d2Likelihoods_
 
VVdouble leafLikelihoods1_
 
VVdouble leafLikelihoods2_
 
double minimumBrLen_
 
Constraint * brLenConstraint_
 
double brLen_
 

Detailed Description

This class is a simplified version of DRHomogeneousTreeLikelihood for 2-Trees.

Definition at line 65 of file DistanceEstimation.h.

Constructor & Destructor Documentation

◆ TwoTreeLikelihood() [1/2]

TwoTreeLikelihood::TwoTreeLikelihood ( const std::string &  seq1,
const std::string &  seq2,
const SiteContainer &  data,
SubstitutionModel model,
DiscreteDistribution *  rDist,
bool  verbose 
)
throw (Exception
)

◆ TwoTreeLikelihood() [2/2]

TwoTreeLikelihood::TwoTreeLikelihood ( const TwoTreeLikelihood lik)

Definition at line 120 of file DistanceEstimation.cpp.

◆ ~TwoTreeLikelihood()

TwoTreeLikelihood::~TwoTreeLikelihood ( )
virtual

Definition at line 180 of file DistanceEstimation.cpp.

References brLenConstraint_, and shrunkData_.

Member Function Documentation

◆ applyParameters()

void TwoTreeLikelihood::applyParameters ( )
throw (Exception
)
protectedvirtual

All parameters are stores in a parameter list.

This function apply these parameters to the substitution model, to the rate distribution and to the branch lengths.

Definition at line 297 of file DistanceEstimation.cpp.

References brLen_, model_, and bpp::AbstractDiscreteRatesAcrossSitesTreeLikelihood::rateDistribution_.

Referenced by fireParameterChanged().

◆ clone()

TwoTreeLikelihood* bpp::TwoTreeLikelihood::clone ( ) const
inlinevirtual

Implements bpp::TreeLikelihood.

Definition at line 127 of file DistanceEstimation.h.

References TwoTreeLikelihood().

◆ computeTreeD2Likelihood()

◆ computeTreeDLikelihood()

◆ computeTreeLikelihood()

◆ displayLikelihoodArray()

void AbstractDiscreteRatesAcrossSitesTreeLikelihood::displayLikelihoodArray ( const VVVdouble &  likelihoodArray)
staticinherited

◆ enableDerivatives()

void bpp::AbstractTreeLikelihood::enableDerivatives ( bool  yn)
inlinevirtualinherited

Tell if derivatives must be computed.

This methods calls the enableFirstOrderDerivatives and enableSecondOrderDerivatives.

Parameters
ynYes or no.

Implements bpp::TreeLikelihood.

Definition at line 292 of file AbstractTreeLikelihood.h.

References bpp::AbstractTreeLikelihood::computeFirstOrderDerivatives_, and bpp::AbstractTreeLikelihood::computeSecondOrderDerivatives_.

Referenced by bpp::AbstractDiscreteRatesAcrossSitesTreeLikelihood::AbstractDiscreteRatesAcrossSitesTreeLikelihood(), and bpp::DistanceEstimation::computeMatrix().

◆ enableFirstOrderDerivatives() [1/2]

void bpp::AbstractTreeLikelihood::enableFirstOrderDerivatives ( bool  yn)
inlineinherited

◆ enableFirstOrderDerivatives() [2/2]

bool bpp::AbstractTreeLikelihood::enableFirstOrderDerivatives ( ) const
inlineinherited

◆ enableSecondOrderDerivatives() [1/2]

void bpp::AbstractTreeLikelihood::enableSecondOrderDerivatives ( bool  yn)
inlineinherited

◆ enableSecondOrderDerivatives() [2/2]

bool bpp::AbstractTreeLikelihood::enableSecondOrderDerivatives ( ) const
inlineinherited

◆ fireParameterChanged()

◆ getAlphabet()

const Alphabet* bpp::AbstractTreeLikelihood::getAlphabet ( ) const
inlinevirtualinherited

Get the alphabet associated to the dataset.

Returns
the alphabet associated to the dataset.

Implements bpp::TreeLikelihood.

Definition at line 285 of file AbstractTreeLikelihood.h.

References bpp::AbstractTreeLikelihood::data_.

◆ getAlphabetStateAsChar()

std::string bpp::TwoTreeLikelihood::getAlphabetStateAsChar ( size_t  i) const
inlinevirtual
Returns
the alphabet state corresponding to the given model state.

Implements bpp::TreeLikelihood.

Definition at line 146 of file DistanceEstimation.h.

References bpp::SubstitutionModel::getAlphabetStateAsChar(), and model_.

◆ getAlphabetStateAsInt()

int bpp::TwoTreeLikelihood::getAlphabetStateAsInt ( size_t  i) const
inlinevirtual
Returns
the alphabet state corresponding to the given model state.

Implements bpp::TreeLikelihood.

Definition at line 144 of file DistanceEstimation.h.

References bpp::SubstitutionModel::getAlphabetStateAsInt(), and model_.

◆ getAlphabetStates()

const std::vector<int>& bpp::TwoTreeLikelihood::getAlphabetStates ( ) const
inlinevirtual
Returns
the alphabet states corresponding to all model states.

Implements bpp::TreeLikelihood.

Definition at line 142 of file DistanceEstimation.h.

References bpp::SubstitutionModel::getAlphabetStates(), and model_.

◆ getBranchLengthsParameters()

ParameterList TwoTreeLikelihood::getBranchLengthsParameters ( ) const
virtual

Get the branch lengths parameters.

Returns
A ParameterList with all branch lengths.

Implements bpp::TreeLikelihood.

Definition at line 197 of file DistanceEstimation.cpp.

References brLenParameters_, and bpp::AbstractTreeLikelihood::initialized_.

Referenced by bpp::DistanceEstimation::computeMatrix().

◆ getData()

const SiteContainer* bpp::AbstractTreeLikelihood::getData ( ) const
inlinevirtualinherited

Get the dataset for which the likelihood must be evaluated.

Returns
A pointer toward the site container where the sequences are stored.

Implements bpp::TreeLikelihood.

Definition at line 284 of file AbstractTreeLikelihood.h.

References bpp::AbstractTreeLikelihood::data_.

Referenced by bpp::RNonHomogeneousMixedTreeLikelihood::init().

◆ getDerivableParameters()

ParameterList AbstractDiscreteRatesAcrossSitesTreeLikelihood::getDerivableParameters ( ) const
virtualinherited

All derivable parameters.

Usually, this contains all branch lengths parameters.

Returns
A ParameterList.

Implements bpp::TreeLikelihood.

Reimplemented in bpp::RHomogeneousClockTreeLikelihood.

Definition at line 73 of file AbstractDiscreteRatesAcrossSitesTreeLikelihood.cpp.

◆ getFirstOrderDerivative()

double TwoTreeLikelihood::getFirstOrderDerivative ( const std::string &  variable) const
throw (Exception
)

Definition at line 590 of file DistanceEstimation.cpp.

◆ getLikelihood()

double TwoTreeLikelihood::getLikelihood ( ) const
virtual

Get the likelihood for the whole dataset.

Returns
The likelihood of the dataset.

Implements bpp::TreeLikelihood.

Definition at line 213 of file DistanceEstimation.cpp.

References nbDistinctSites_, rootLikelihoodsSR_, and rootWeights_.

◆ getLikelihoodData() [1/2]

TreeLikelihoodData* bpp::TwoTreeLikelihood::getLikelihoodData ( )
throw (NotImplementedException
)
inlinevirtual
Returns
The underlying likelihood data structure.

Implements bpp::TreeLikelihood.

Definition at line 148 of file DistanceEstimation.h.

◆ getLikelihoodData() [2/2]

const TreeLikelihoodData* bpp::TwoTreeLikelihood::getLikelihoodData ( ) const
throw (NotImplementedException
)
inlinevirtual
Returns
The underlying likelihood data structure.

Implements bpp::TreeLikelihood.

Definition at line 152 of file DistanceEstimation.h.

◆ getLikelihoodForASite()

double TwoTreeLikelihood::getLikelihoodForASite ( size_t  site) const
virtual

Get the likelihood for a site.

Parameters
siteThe site index to analyse.
Returns
The likelihood for site site.

Implements bpp::TreeLikelihood.

Definition at line 237 of file DistanceEstimation.cpp.

References rootLikelihoodsSR_, and rootPatternLinks_.

◆ getLikelihoodForASiteForARateClass()

double TwoTreeLikelihood::getLikelihoodForASiteForARateClass ( size_t  site,
size_t  rateClass 
) const
virtual

Get the likelihood for a site knowing its rate class.

Parameters
siteThe site index.
rateClassThe rate class index.
Returns
The likelihood for the specified site and rate class.

Implements bpp::DiscreteRatesAcrossSitesTreeLikelihood.

Definition at line 251 of file DistanceEstimation.cpp.

References rootLikelihoodsS_, and rootPatternLinks_.

◆ getLikelihoodForASiteForARateClassForAState()

double TwoTreeLikelihood::getLikelihoodForASiteForARateClassForAState ( size_t  site,
size_t  rateClass,
int  state 
) const
virtual

Get the likelihood for a site knowing its rate class and its ancestral state.

Parameters
siteThe site index.
rateClassThe rate class index.
stateThe ancestral state.
Returns
The likelihood for the specified site and rate class and ancestral state.

Implements bpp::DiscreteRatesAcrossSitesTreeLikelihood.

Definition at line 265 of file DistanceEstimation.cpp.

References rootLikelihoods_, and rootPatternLinks_.

◆ getLikelihoodForASiteForAState()

double AbstractDiscreteRatesAcrossSitesTreeLikelihood::getLikelihoodForASiteForAState ( size_t  site,
int  state 
) const
virtualinherited

Get the likelihood for a site and for a state.

Parameters
siteThe site index to analyse.
stateThe state to consider.
Returns
The likelihood for site site and state state.

Implements bpp::TreeLikelihood.

Definition at line 111 of file AbstractDiscreteRatesAcrossSitesTreeLikelihood.cpp.

◆ getLikelihoodForEachSite()

Vdouble AbstractTreeLikelihood::getLikelihoodForEachSite ( ) const
virtualinherited

Get the likelihood for each site.

Returns
A vector with all likelihoods for each site.

Implements bpp::TreeLikelihood.

Definition at line 46 of file AbstractTreeLikelihood.cpp.

References bpp::TreeLikelihood::getLikelihoodForASite(), and bpp::AbstractTreeLikelihood::getNumberOfSites().

◆ getLikelihoodForEachSiteForEachRateClass()

VVdouble AbstractDiscreteRatesAcrossSitesTreeLikelihood::getLikelihoodForEachSiteForEachRateClass ( ) const
virtualinherited

Get the likelihood for each site and each rate class.

Returns
A two-dimension vector with all likelihoods.

Implements bpp::DiscreteRatesAcrossSitesTreeLikelihood.

Definition at line 93 of file AbstractDiscreteRatesAcrossSitesTreeLikelihood.cpp.

◆ getLikelihoodForEachSiteForEachRateClassForEachState()

VVVdouble AbstractDiscreteRatesAcrossSitesTreeLikelihood::getLikelihoodForEachSiteForEachRateClassForEachState ( ) const
virtualinherited

Get the likelihood for each site and each rate class and each state.

Returns
A three-dimension vector with all likelihoods.

Implements bpp::DiscreteRatesAcrossSitesTreeLikelihood.

Definition at line 156 of file AbstractDiscreteRatesAcrossSitesTreeLikelihood.cpp.

◆ getLikelihoodForEachSiteForEachState()

VVdouble AbstractTreeLikelihood::getLikelihoodForEachSiteForEachState ( ) const
virtualinherited

Get the likelihood for each site and for each state.

Returns
A 2d vector with all likelihoods for each site and for each state.

Implements bpp::TreeLikelihood.

Definition at line 64 of file AbstractTreeLikelihood.cpp.

References bpp::TreeLikelihood::getLikelihoodForASiteForAState(), bpp::AbstractTreeLikelihood::getNumberOfSites(), and bpp::TreeLikelihood::getNumberOfStates().

◆ getLogLikelihood()

double TwoTreeLikelihood::getLogLikelihood ( ) const
virtual

Get the logarithm of the likelihood for the whole dataset.

Returns
The logarithm of the likelihood of the dataset.

Implements bpp::TreeLikelihood.

Definition at line 225 of file DistanceEstimation.cpp.

References nbDistinctSites_, rootLikelihoodsSR_, and rootWeights_.

Referenced by getValue().

◆ getLogLikelihoodForASite()

double TwoTreeLikelihood::getLogLikelihoodForASite ( size_t  site) const
virtual

Get the logarithm of the likelihood for a site.

Parameters
siteThe site index to analyse.
Returns
The logarithm of the likelihood for site site.

Implements bpp::TreeLikelihood.

Definition at line 244 of file DistanceEstimation.cpp.

References rootLikelihoodsSR_, and rootPatternLinks_.

◆ getLogLikelihoodForASiteForARateClass()

double TwoTreeLikelihood::getLogLikelihoodForASiteForARateClass ( size_t  site,
size_t  rateClass 
) const
virtual

Get the logarithm of the likelihood for a site knowing its rate class.

Parameters
siteThe site index.
rateClassThe rate class index.
Returns
The logarithm of the likelihood for the specified site and rate class.

Implements bpp::DiscreteRatesAcrossSitesTreeLikelihood.

Definition at line 258 of file DistanceEstimation.cpp.

References rootLikelihoodsS_, and rootPatternLinks_.

◆ getLogLikelihoodForASiteForARateClassForAState()

double TwoTreeLikelihood::getLogLikelihoodForASiteForARateClassForAState ( size_t  site,
size_t  rateClass,
int  state 
) const
virtual

Get the logarithm of the likelihood for a site knowing its rate class and its ancestral state.

Parameters
siteThe site index.
rateClassThe rate class index.
stateThe ancestral state.
Returns
The logarithm of the likelihood for the specified site and rate class and ancestral state..

Implements bpp::DiscreteRatesAcrossSitesTreeLikelihood.

Definition at line 272 of file DistanceEstimation.cpp.

References rootLikelihoods_, and rootPatternLinks_.

◆ getLogLikelihoodForASiteForAState()

double AbstractDiscreteRatesAcrossSitesTreeLikelihood::getLogLikelihoodForASiteForAState ( size_t  site,
int  state 
) const
virtualinherited

Get the logarithm of the likelihood for a site and for a state.

Parameters
siteThe site index to analyse.
stateThe state to consider.
Returns
The logarithm of the likelihood for site site and state state.

Implements bpp::TreeLikelihood.

Definition at line 124 of file AbstractDiscreteRatesAcrossSitesTreeLikelihood.cpp.

◆ getLogLikelihoodForEachSite()

Vdouble AbstractTreeLikelihood::getLogLikelihoodForEachSite ( ) const
virtualinherited

Get the logarithm of the likelihood for each site.

Returns
A vector with all log likelihoods for each site.

Implements bpp::TreeLikelihood.

Definition at line 55 of file AbstractTreeLikelihood.cpp.

References bpp::TreeLikelihood::getLogLikelihoodForASite(), and bpp::AbstractTreeLikelihood::getNumberOfSites().

◆ getLogLikelihoodForEachSiteForEachRateClass()

VVdouble AbstractDiscreteRatesAcrossSitesTreeLikelihood::getLogLikelihoodForEachSiteForEachRateClass ( ) const
virtualinherited

Get the logarithm of the likelihood for each site and each rate class.

Returns
A two-dimension vector with all log likelihoods: V[i][j] = likelihood of site i and rate class j.

Implements bpp::DiscreteRatesAcrossSitesTreeLikelihood.

Definition at line 138 of file AbstractDiscreteRatesAcrossSitesTreeLikelihood.cpp.

◆ getLogLikelihoodForEachSiteForEachRateClassForEachState()

VVVdouble AbstractDiscreteRatesAcrossSitesTreeLikelihood::getLogLikelihoodForEachSiteForEachRateClassForEachState ( ) const
virtualinherited

Get the logarithm of the likelihood for each site and each rate class and each state.

Returns
A three-dimension vector with all log likelihoods: V[i][j][k} = likelihood of site i and rate class j and state k.

Implements bpp::DiscreteRatesAcrossSitesTreeLikelihood.

Definition at line 179 of file AbstractDiscreteRatesAcrossSitesTreeLikelihood.cpp.

◆ getLogLikelihoodForEachSiteForEachState()

VVdouble AbstractTreeLikelihood::getLogLikelihoodForEachSiteForEachState ( ) const
virtualinherited

Get the logarithm of the likelihood for each site and for each state.

Returns
A 2d vector with all log likelihoods for each site and for each state.

Implements bpp::TreeLikelihood.

Definition at line 81 of file AbstractTreeLikelihood.cpp.

References bpp::TreeLikelihood::getLogLikelihoodForASiteForAState(), bpp::AbstractTreeLikelihood::getNumberOfSites(), and bpp::TreeLikelihood::getNumberOfStates().

◆ getMinimumBranchLength()

virtual double bpp::TwoTreeLikelihood::getMinimumBranchLength ( ) const
inlinevirtual

Definition at line 254 of file DistanceEstimation.h.

References minimumBrLen_.

Referenced by bpp::DistanceEstimation::computeMatrix().

◆ getNewBranchModelIterator()

ConstBranchModelIterator* bpp::TwoTreeLikelihood::getNewBranchModelIterator ( int  nodeId) const
throw (NotImplementedException
)
inlinevirtual

Implements bpp::TreeLikelihood.

Definition at line 199 of file DistanceEstimation.h.

◆ getNewSiteModelIterator()

ConstSiteModelIterator* bpp::TwoTreeLikelihood::getNewSiteModelIterator ( size_t  siteIndex) const
throw (NotImplementedException
)
inlinevirtual

Implements bpp::TreeLikelihood.

Definition at line 204 of file DistanceEstimation.h.

◆ getNonDerivableParameters()

ParameterList AbstractDiscreteRatesAcrossSitesTreeLikelihood::getNonDerivableParameters ( ) const
virtualinherited

All non derivable parameters.

Usually, this contains all substitution model parameters and rate distribution.

Returns
A ParameterList.

Implements bpp::TreeLikelihood.

Reimplemented in bpp::RHomogeneousClockTreeLikelihood.

Definition at line 82 of file AbstractDiscreteRatesAcrossSitesTreeLikelihood.cpp.

◆ getNumberOfClasses()

size_t bpp::AbstractDiscreteRatesAcrossSitesTreeLikelihood::getNumberOfClasses ( ) const
inlinevirtualinherited

Get the number of classes.

Returns
The Number of classes.

Implements bpp::DiscreteRatesAcrossSitesTreeLikelihood.

Definition at line 108 of file AbstractDiscreteRatesAcrossSitesTreeLikelihood.h.

References bpp::AbstractDiscreteRatesAcrossSitesTreeLikelihood::rateDistribution_.

◆ getNumberOfSites()

size_t bpp::AbstractTreeLikelihood::getNumberOfSites ( ) const
inlinevirtualinherited

◆ getNumberOfStates()

size_t bpp::TwoTreeLikelihood::getNumberOfStates ( ) const
inlinevirtual
Returns
the number of model states of the underlying Markov chain.

Implements bpp::TreeLikelihood.

Definition at line 140 of file DistanceEstimation.h.

References bpp::SubstitutionModel::getNumberOfStates(), and model_.

◆ getPosteriorProbabilitiesOfEachRate()

VVdouble AbstractDiscreteRatesAcrossSitesTreeLikelihood::getPosteriorProbabilitiesOfEachRate ( ) const
virtualinherited

Get the posterior probability for each site of belonging to a particular rate class.

Returns
A two-dimension vector with all posterior probabilities: V[i][j] = probablity for site i of belonging to rate class j.

Implements bpp::DiscreteRatesAcrossSitesTreeLikelihood.

Definition at line 202 of file AbstractDiscreteRatesAcrossSitesTreeLikelihood.cpp.

◆ getPosteriorRateOfEachSite()

Vdouble AbstractDiscreteRatesAcrossSitesTreeLikelihood::getPosteriorRateOfEachSite ( ) const
virtualinherited

Get the posterior rate, i.e. averaged over all classes and weighted with posterior probabilities, for each site.

Returns
A vector with all rates.

Implements bpp::DiscreteRatesAcrossSitesTreeLikelihood.

Definition at line 220 of file AbstractDiscreteRatesAcrossSitesTreeLikelihood.cpp.

◆ getRateClassWithMaxPostProbOfEachSite()

vector< size_t > AbstractDiscreteRatesAcrossSitesTreeLikelihood::getRateClassWithMaxPostProbOfEachSite ( ) const
virtualinherited

Get the posterior rate class (the one with maximum posterior probability) for each site.

Returns
A vector with all rate classes indexes.

Implements bpp::DiscreteRatesAcrossSitesTreeLikelihood.

Definition at line 239 of file AbstractDiscreteRatesAcrossSitesTreeLikelihood.cpp.

◆ getRateDistribution() [1/2]

const DiscreteDistribution* bpp::AbstractDiscreteRatesAcrossSitesTreeLikelihood::getRateDistribution ( ) const
inlinevirtualinherited

Get the rate distribution used for the computation.

Returns
A const pointer toward the rate distribution of this instance.

Implements bpp::DiscreteRatesAcrossSitesTreeLikelihood.

Definition at line 106 of file AbstractDiscreteRatesAcrossSitesTreeLikelihood.h.

References bpp::AbstractDiscreteRatesAcrossSitesTreeLikelihood::rateDistribution_.

◆ getRateDistribution() [2/2]

DiscreteDistribution* bpp::AbstractDiscreteRatesAcrossSitesTreeLikelihood::getRateDistribution ( )
inlinevirtualinherited

Get the rate distribution used for the computation.

Returns
A pointer toward the rate distribution of this instance.

Implements bpp::DiscreteRatesAcrossSitesTreeLikelihood.

Definition at line 107 of file AbstractDiscreteRatesAcrossSitesTreeLikelihood.h.

References bpp::AbstractDiscreteRatesAcrossSitesTreeLikelihood::rateDistribution_.

◆ getRateDistributionParameters()

ParameterList AbstractDiscreteRatesAcrossSitesTreeLikelihood::getRateDistributionParameters ( ) const
virtualinherited

Get the parameters associated to the rate distirbution.

Returns
A ParameterList object with all rate distribution parameters.

Implements bpp::DiscreteRatesAcrossSitesTreeLikelihood.

Reimplemented in bpp::AbstractNonHomogeneousTreeLikelihood, and bpp::AbstractHomogeneousTreeLikelihood.

Definition at line 64 of file AbstractDiscreteRatesAcrossSitesTreeLikelihood.cpp.

Referenced by bpp::AbstractHomogeneousTreeLikelihood::getRateDistributionParameters(), and bpp::AbstractNonHomogeneousTreeLikelihood::getRateDistributionParameters().

◆ getRateWithMaxPostProbOfEachSite()

Vdouble AbstractDiscreteRatesAcrossSitesTreeLikelihood::getRateWithMaxPostProbOfEachSite ( ) const
virtualinherited

Get the posterior rate (the one with maximum posterior probability) for each site.

Returns
A vector with all rate classes indexes.

Implements bpp::DiscreteRatesAcrossSitesTreeLikelihood.

Definition at line 253 of file AbstractDiscreteRatesAcrossSitesTreeLikelihood.cpp.

◆ getRootFrequencies()

const std::vector<double>& bpp::TwoTreeLikelihood::getRootFrequencies ( size_t  siteIndex) const
inlinevirtual

Get the values of the frequencies for each state in the alphabet at the root node.

For reversible models, these are the equilibrium frequencies. For non-reversible models, these usually are distinct parameters.

For models without site partitioning, the set of frequencies is the same for all positions. For partition models, the frequencies may differ from one site to another.

Parameters
siteIndexThe index of the alignment position.
See also
getSiteIndex
Returns
A vector with ancestral frequencies for each state in the alphabet;

Implements bpp::TreeLikelihood.

Definition at line 164 of file DistanceEstimation.h.

References bpp::SubstitutionModel::getFrequencies(), and model_.

◆ getSecondOrderDerivative() [1/2]

double TwoTreeLikelihood::getSecondOrderDerivative ( const std::string &  variable) const
throw (Exception
)

Definition at line 621 of file DistanceEstimation.cpp.

◆ getSecondOrderDerivative() [2/2]

double bpp::TwoTreeLikelihood::getSecondOrderDerivative ( const std::string &  variable1,
const std::string &  variable2 
) const
throw (Exception
)
inline

Definition at line 241 of file DistanceEstimation.h.

◆ getSiteIndex()

size_t bpp::TwoTreeLikelihood::getSiteIndex ( size_t  site) const
throw (IndexOutOfBoundsException
)
inlinevirtual

Get the index (used for inner computations) of a given site (original alignment column).

Parameters
siteAn alignment position.
Returns
The site index corresponding to the given input alignment position.

Implements bpp::TreeLikelihood.

Definition at line 165 of file DistanceEstimation.h.

References rootPatternLinks_.

◆ getSubstitutionModel() [1/4]

SubstitutionModel* bpp::TwoTreeLikelihood::getSubstitutionModel ( int  nodeId,
size_t  siteIndex 
)
throw (NodeNotFoundException
)
inlinevirtual

Get the substitution model associated to a given node and alignment column.

Parameters
nodeIdThe id of the request node.
siteIndexThe index of the alignment position.
See also
getSiteIndex
Returns
A pointer toward the corresponding model.
Exceptions
NodeNotFoundExceptionThis exception may be thrown if the node is not found (depending on the implementation).

Implements bpp::TreeLikelihood.

Definition at line 162 of file DistanceEstimation.h.

References model_.

◆ getSubstitutionModel() [2/4]

const SubstitutionModel* bpp::TwoTreeLikelihood::getSubstitutionModel ( int  nodeId,
size_t  siteIndex 
) const
throw (NodeNotFoundException
)
inlinevirtual

Get the substitution model associated to a given node and alignment column.

Parameters
nodeIdThe id of the request node.
siteIndexThe index of the alignment position.
See also
getSiteIndex
Returns
A pointer toward the corresponding model.
Exceptions
NodeNotFoundExceptionThis exception may be thrown if the node is not found (depending on the implementation).

Implements bpp::TreeLikelihood.

Definition at line 163 of file DistanceEstimation.h.

References model_.

◆ getSubstitutionModel() [3/4]

const SubstitutionModel* bpp::TwoTreeLikelihood::getSubstitutionModel ( ) const
inline

Get the substitution model used for the computation.

Returns
A const pointer toward the substitution model of this instance.

Definition at line 190 of file DistanceEstimation.h.

References model_.

◆ getSubstitutionModel() [4/4]

SubstitutionModel* bpp::TwoTreeLikelihood::getSubstitutionModel ( )
inline

Get the substitution model used for the computation.

Returns
A pointer toward the substitution model of this instance.

Definition at line 197 of file DistanceEstimation.h.

References model_.

◆ getSubstitutionModelParameters()

ParameterList TwoTreeLikelihood::getSubstitutionModelParameters ( ) const
virtual

Get the parameters associated to substitution model(s).

Returns
A ParameterList.

Implements bpp::TreeLikelihood.

Definition at line 205 of file DistanceEstimation.cpp.

References bpp::AbstractTreeLikelihood::initialized_, and model_.

◆ getTransitionProbabilities()

VVdouble AbstractDiscreteRatesAcrossSitesTreeLikelihood::getTransitionProbabilities ( int  nodeId,
size_t  siteIndex 
) const
virtualinherited

Retrieves all Pij(t) for a particular branch, defined by the upper node and site.

These intermediate results may be used by other methods.

Parameters
nodeIdThe node defining the branch of interest.
siteIndexThe index of the alignment position.
See also
getSiteIndex
Returns
An array of dimension 2, where a[x][y] is the probability of substituting from x to y.

Implements bpp::TreeLikelihood.

Definition at line 311 of file AbstractDiscreteRatesAcrossSitesTreeLikelihood.cpp.

◆ getTransitionProbabilitiesPerRateClass()

VVVdouble bpp::TwoTreeLikelihood::getTransitionProbabilitiesPerRateClass ( int  nodeId,
size_t  siteIndex 
) const
inlinevirtual

This method is not applicable for this object.

Implements bpp::DiscreteRatesAcrossSitesTreeLikelihood.

Definition at line 169 of file DistanceEstimation.h.

References pxy_.

◆ getTree()

const Tree& bpp::AbstractTreeLikelihood::getTree ( ) const
inlinevirtualinherited

Get the tree (topology and branch lengths).

Returns
The tree of this TreeLikelihood object.

Implements bpp::TreeLikelihood.

Definition at line 291 of file AbstractTreeLikelihood.h.

References bpp::AbstractTreeLikelihood::tree_.

Referenced by bpp::NNIHomogeneousTreeLikelihood::getTopology(), and bpp::RNonHomogeneousMixedTreeLikelihood::init().

◆ getValue()

double TwoTreeLikelihood::getValue ( ) const
throw (Exception
)

Definition at line 411 of file DistanceEstimation.cpp.

References getLogLikelihood().

◆ initBranchLengthsParameters()

void TwoTreeLikelihood::initBranchLengthsParameters ( )
virtual

◆ initialize()

void TwoTreeLikelihood::initialize ( )
throw (Exception
)
virtual

Init the likelihood object.

This method is used to initialize all parameters. It is typically called after the constructor and the setData method. It contains virtual methods that can't be called in the constructor.

Exceptions
Exceptionif something bad happened, for instance if no data are associated to the likelihood function.

Implements bpp::TreeLikelihood.

Definition at line 188 of file DistanceEstimation.cpp.

References fireParameterChanged(), bpp::AbstractTreeLikelihood::initialized_, and initParameters().

Referenced by bpp::DistanceEstimation::computeMatrix().

◆ initParameters()

void TwoTreeLikelihood::initParameters ( )
protectedvirtual

This builds the parameters list from all parametrizable objects, i.e. substitution model, rate distribution and tree.

Definition at line 279 of file DistanceEstimation.cpp.

References brLenParameters_, initBranchLengthsParameters(), model_, and bpp::AbstractDiscreteRatesAcrossSitesTreeLikelihood::rateDistribution_.

Referenced by initialize().

◆ initTreeLikelihoods()

void TwoTreeLikelihood::initTreeLikelihoods ( const SequenceContainer &  sequences)
throw (Exception
)
protectedvirtual

This method initializes the leaves according to a sequence container.

Here the container shrunkData_ is used. Likelihood is set to 1 for the state corresponding to the sequence site, otherwise it is set to 0.

The two likelihood arrays are initialized according to alphabet size and sequences length, and filled with 1.

NB: This method is recursive.

Parameters
sequencesThe sequence container to use.

Definition at line 419 of file DistanceEstimation.cpp.

◆ isInitialized()

bool bpp::AbstractTreeLikelihood::isInitialized ( ) const
inlinevirtualinherited

◆ operator=()

◆ resetLikelihoodArray()

◆ setData()

void bpp::TwoTreeLikelihood::setData ( const SiteContainer &  sites)
throw (Exception
)
inlinevirtual

Set the dataset for which the likelihood must be evaluated.

Parameters
sitesThe data set to use.

Implements bpp::TreeLikelihood.

Definition at line 170 of file DistanceEstimation.h.

◆ setMinimumBranchLength()

virtual void bpp::TwoTreeLikelihood::setMinimumBranchLength ( double  minimum)
inlinevirtual

◆ setParameters()

void TwoTreeLikelihood::setParameters ( const ParameterList &  parameters)
throw (ParameterNotFoundException,
ConstraintException
)

Implements the Function interface.

Update the parameter list and call the applyParameters() method. Then compute the likelihoods at each node (computeLikelihood() method) and call the getLogLikelihood() method.

If a subset of the whole parameter list is passed to the function, only these parameters are updated and the other remain constant (i.e. equal to their last value).

Parameters
parametersThe parameter list to pass to the function.

Definition at line 322 of file DistanceEstimation.cpp.

Member Data Documentation

◆ brLen_

double bpp::TwoTreeLikelihood::brLen_
private

◆ brLenConstraint_

Constraint* bpp::TwoTreeLikelihood::brLenConstraint_
private

◆ brLenParameters_

ParameterList bpp::TwoTreeLikelihood::brLenParameters_
private

◆ computeFirstOrderDerivatives_

◆ computeSecondOrderDerivatives_

◆ d2Likelihoods_

Vdouble bpp::TwoTreeLikelihood::d2Likelihoods_
mutableprivate

Definition at line 108 of file DistanceEstimation.h.

Referenced by computeTreeD2Likelihood(), and operator=().

◆ d2pxy_

VVVdouble bpp::TwoTreeLikelihood::d2pxy_
mutableprivate

Definition at line 78 of file DistanceEstimation.h.

Referenced by computeTreeD2Likelihood(), fireParameterChanged(), and operator=().

◆ data_

◆ dLikelihoods_

Vdouble bpp::TwoTreeLikelihood::dLikelihoods_
mutableprivate

Definition at line 107 of file DistanceEstimation.h.

Referenced by computeTreeDLikelihood(), and operator=().

◆ dpxy_

VVVdouble bpp::TwoTreeLikelihood::dpxy_
mutableprivate

Definition at line 76 of file DistanceEstimation.h.

Referenced by computeTreeDLikelihood(), fireParameterChanged(), and operator=().

◆ initialized_

◆ leafLikelihoods1_

VVdouble bpp::TwoTreeLikelihood::leafLikelihoods1_
mutableprivate

◆ leafLikelihoods2_

VVdouble bpp::TwoTreeLikelihood::leafLikelihoods2_
mutableprivate

◆ minimumBrLen_

double bpp::TwoTreeLikelihood::minimumBrLen_
private

◆ model_

◆ nbClasses_

size_t bpp::TwoTreeLikelihood::nbClasses_
private

◆ nbDistinctSites_

size_t bpp::TwoTreeLikelihood::nbDistinctSites_
private

◆ nbSites_

size_t bpp::TwoTreeLikelihood::nbSites_
private

Definition at line 99 of file DistanceEstimation.h.

Referenced by operator=().

◆ nbStates_

size_t bpp::TwoTreeLikelihood::nbStates_
private

◆ pxy_

VVVdouble bpp::TwoTreeLikelihood::pxy_
mutableprivate

◆ rateDistribution_

DiscreteDistribution* bpp::AbstractDiscreteRatesAcrossSitesTreeLikelihood::rateDistribution_
protectedinherited

Definition at line 61 of file AbstractDiscreteRatesAcrossSitesTreeLikelihood.h.

Referenced by bpp::AbstractHomogeneousTreeLikelihood::applyParameters(), bpp::AbstractNonHomogeneousTreeLikelihood::applyParameters(), applyParameters(), bpp::DRHomogeneousTreeLikelihood::computeRootLikelihood(), bpp::DRNonHomogeneousTreeLikelihood::computeRootLikelihood(), bpp::AbstractHomogeneousTreeLikelihood::computeTransitionProbabilitiesForNode(), bpp::RNonHomogeneousMixedTreeLikelihood::computeTransitionProbabilitiesForNode(), bpp::AbstractNonHomogeneousTreeLikelihood::computeTransitionProbabilitiesForNode(), computeTreeD2Likelihood(), bpp::DRHomogeneousTreeLikelihood::computeTreeD2LikelihoodAtNode(), bpp::DRNonHomogeneousTreeLikelihood::computeTreeD2LikelihoodAtNode(), computeTreeDLikelihood(), bpp::DRHomogeneousTreeLikelihood::computeTreeDLikelihoodAtNode(), bpp::DRNonHomogeneousTreeLikelihood::computeTreeDLikelihoodAtNode(), computeTreeLikelihood(), bpp::DRHomogeneousTreeLikelihood::fireParameterChanged(), bpp::DRNonHomogeneousTreeLikelihood::fireParameterChanged(), bpp::RHomogeneousTreeLikelihood::fireParameterChanged(), bpp::RNonHomogeneousTreeLikelihood::fireParameterChanged(), fireParameterChanged(), bpp::RNonHomogeneousMixedTreeLikelihood::fireParameterChanged(), bpp::RNonHomogeneousTreeLikelihood::getD2LikelihoodForASite(), bpp::RHomogeneousTreeLikelihood::getD2LikelihoodForASite(), bpp::RNonHomogeneousTreeLikelihood::getDLikelihoodForASite(), bpp::RHomogeneousTreeLikelihood::getDLikelihoodForASite(), bpp::RHomogeneousTreeLikelihood::getLikelihoodForASite(), bpp::RNonHomogeneousTreeLikelihood::getLikelihoodForASite(), bpp::RHomogeneousTreeLikelihood::getLogLikelihoodForASite(), bpp::RNonHomogeneousTreeLikelihood::getLogLikelihoodForASite(), bpp::AbstractDiscreteRatesAcrossSitesTreeLikelihood::getNumberOfClasses(), bpp::AbstractDiscreteRatesAcrossSitesTreeLikelihood::getRateDistribution(), bpp::RNonHomogeneousMixedTreeLikelihood::init(), bpp::AbstractHomogeneousTreeLikelihood::initParameters(), bpp::AbstractNonHomogeneousTreeLikelihood::initParameters(), initParameters(), and bpp::AbstractDiscreteRatesAcrossSitesTreeLikelihood::operator=().

◆ rootLikelihoods_

VVVdouble bpp::TwoTreeLikelihood::rootLikelihoods_
mutableprivate

◆ rootLikelihoodsS_

VVdouble bpp::TwoTreeLikelihood::rootLikelihoodsS_
mutableprivate

◆ rootLikelihoodsSR_

Vdouble bpp::TwoTreeLikelihood::rootLikelihoodsSR_
mutableprivate

◆ rootPatternLinks_

std::vector<size_t> bpp::TwoTreeLikelihood::rootPatternLinks_
private

As previous, but for the global container.

The size of this vector is equal to the number of sites in the container, each element corresponds to a site in the container and points to the corresponding column in the likelihood array of the root node. If the container contains no repeated site, there will be a strict equivalence between each site and the likelihood array of the root node. However, if this is not the case, some pointers may point toward the same element in the likelihood array.

Definition at line 91 of file DistanceEstimation.h.

Referenced by getLikelihoodForASite(), getLikelihoodForASiteForARateClass(), getLikelihoodForASiteForARateClassForAState(), getLogLikelihoodForASite(), getLogLikelihoodForASiteForARateClass(), getLogLikelihoodForASiteForARateClassForAState(), getSiteIndex(), and operator=().

◆ rootWeights_

std::vector<unsigned int> bpp::TwoTreeLikelihood::rootWeights_
private

The frequency of each site.

Definition at line 96 of file DistanceEstimation.h.

Referenced by getLikelihood(), getLogLikelihood(), and operator=().

◆ seqnames_

std::vector<std::string> bpp::TwoTreeLikelihood::seqnames_
private

Definition at line 70 of file DistanceEstimation.h.

Referenced by operator=().

◆ shrunkData_

SiteContainer* bpp::TwoTreeLikelihood::shrunkData_
private

Definition at line 69 of file DistanceEstimation.h.

Referenced by operator=(), and ~TwoTreeLikelihood().

◆ tree_

TreeTemplate<Node>* bpp::AbstractTreeLikelihood::tree_
mutableprotectedinherited

Definition at line 226 of file AbstractTreeLikelihood.h.

Referenced by bpp::AbstractHomogeneousTreeLikelihood::AbstractHomogeneousTreeLikelihood(), bpp::AbstractNonHomogeneousTreeLikelihood::AbstractNonHomogeneousTreeLikelihood(), bpp::AbstractTreeLikelihood::AbstractTreeLikelihood(), bpp::DRHomogeneousTreeLikelihood::computeLikelihoodAtNode(), bpp::DRNonHomogeneousTreeLikelihood::computeLikelihoodAtNode(), bpp::DRHomogeneousTreeLikelihood::computeRootLikelihood(), bpp::DRNonHomogeneousTreeLikelihood::computeRootLikelihood(), bpp::RNonHomogeneousTreeLikelihood::computeTreeD2Likelihood(), bpp::RNonHomogeneousMixedTreeLikelihood::computeTreeD2Likelihood(), bpp::RNonHomogeneousTreeLikelihood::computeTreeDLikelihood(), bpp::RNonHomogeneousMixedTreeLikelihood::computeTreeDLikelihood(), bpp::DRHomogeneousTreeLikelihood::computeTreeLikelihood(), bpp::DRNonHomogeneousTreeLikelihood::computeTreeLikelihood(), bpp::RNonHomogeneousTreeLikelihood::computeTreeLikelihood(), bpp::RHomogeneousTreeLikelihood::computeTreeLikelihood(), bpp::DRHomogeneousTreeLikelihood::DRHomogeneousTreeLikelihood(), bpp::DRNonHomogeneousTreeLikelihood::DRNonHomogeneousTreeLikelihood(), bpp::DRNonHomogeneousTreeLikelihood::fireParameterChanged(), bpp::RNonHomogeneousTreeLikelihood::fireParameterChanged(), bpp::RNonHomogeneousMixedTreeLikelihood::fireParameterChanged(), bpp::RNonHomogeneousTreeLikelihood::getD2LikelihoodForASiteForARateClass(), bpp::RHomogeneousTreeLikelihood::getD2LikelihoodForASiteForARateClass(), bpp::RNonHomogeneousTreeLikelihood::getDLikelihoodForASiteForARateClass(), bpp::RHomogeneousTreeLikelihood::getDLikelihoodForASiteForARateClass(), bpp::RHomogeneousTreeLikelihood::getLikelihoodForASiteForARateClass(), bpp::RNonHomogeneousTreeLikelihood::getLikelihoodForASiteForARateClass(), bpp::RHomogeneousTreeLikelihood::getLikelihoodForASiteForARateClassForAState(), bpp::RNonHomogeneousTreeLikelihood::getLikelihoodForASiteForARateClassForAState(), bpp::RHomogeneousTreeLikelihood::getLogLikelihoodForASiteForARateClass(), bpp::RNonHomogeneousTreeLikelihood::getLogLikelihoodForASiteForARateClass(), bpp::RHomogeneousTreeLikelihood::getLogLikelihoodForASiteForARateClassForAState(), bpp::RNonHomogeneousTreeLikelihood::getLogLikelihoodForASiteForARateClassForAState(), bpp::AbstractHomogeneousTreeLikelihood::getNewSiteModelIterator(), bpp::AbstractTreeLikelihood::getTree(), bpp::DRHomogeneousTreeLikelihood::operator=(), bpp::DRNonHomogeneousTreeLikelihood::operator=(), bpp::AbstractHomogeneousTreeLikelihood::operator=(), bpp::RHomogeneousTreeLikelihood::operator=(), bpp::RNonHomogeneousTreeLikelihood::operator=(), bpp::AbstractNonHomogeneousTreeLikelihood::operator=(), bpp::AbstractTreeLikelihood::operator=(), bpp::RHomogeneousTreeLikelihood::RHomogeneousTreeLikelihood(), bpp::RNonHomogeneousTreeLikelihood::RNonHomogeneousTreeLikelihood(), and bpp::AbstractTreeLikelihood::~AbstractTreeLikelihood().


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