46 #include "../PatternTools.h" 48 #include <Bpp/Numeric/VectorTools.h> 49 #include <Bpp/App/ApplicationTools.h> 57 DiscreteDistribution* rDist,
60 bool rootArray)
throw (Exception) :
62 treeLikelihoodsContainer_(),
68 if ((mixedmodel = dynamic_cast<MixedSubstitutionModel*>(model_)) == NULL)
69 throw Exception(
"Bad model: DRHomogeneousMixedTreeLikelihood needs a MixedSubstitutionModel.");
72 for (
size_t i = 0; i < s; i++)
74 treeLikelihoodsContainer_.push_back(
82 const SiteContainer& data,
84 DiscreteDistribution* rDist,
90 treeLikelihoodsContainer_(),
96 if ((mixedmodel = dynamic_cast<MixedSubstitutionModel*>(model_)) == NULL)
97 throw Exception(
"Bad model: DRHomogeneousMixedTreeLikelihood needs a MixedSubstitutionModel.");
101 for (
size_t i = 0; i < s; i++)
103 treeLikelihoodsContainer_.push_back(
115 treeLikelihoodsContainer_.clear();
118 for (
unsigned int i = 0; i < treeLikelihoodsContainer_.size(); i++)
121 probas_.push_back(lik.
probas_[i]);
131 treeLikelihoodsContainer_(lik.treeLikelihoodsContainer_.size()),
132 probas_(lik.probas_.size()),
133 rootArray_(lik.rootArray_)
166 for (
unsigned int i = 0; i < treeLikelihoodsContainer_.size(); i++)
168 treeLikelihoodsContainer_[i]->setData(sites);
182 for (
size_t i = 0; i < s; i++)
186 pl.addParameters(pm->getParameters());
187 pl.includeParameters(getParameters());
220 vector<Vdouble*> llik;
233 x += (*llik[j])[i] *
probas_[j];
235 l *= std::pow(x, (
int)(*w)[i]);
244 vector<Vdouble*> llik;
258 x += (*llik[j])[i] *
probas_[j];
260 la[i] = (*w)[i] * log(x);
262 sort(la.begin(), la.end());
356 VVdouble* likelihoodArray_i = &likelihoodArray[i];
359 Vdouble* likelihoodArray_i_c = &(*likelihoodArray_i)[c];
362 (*likelihoodArray_i_c)[x] = 0;
373 VVdouble* likelihoodArray_i = &likelihoodArray[i];
374 VVdouble* lArray_i = &lArray[i];
378 Vdouble* likelihoodArray_i_c = &(*likelihoodArray_i)[c];
379 Vdouble* lArray_i_c = &(*lArray_i)[c];
381 (*likelihoodArray_i_c)[x] += (*lArray_i_c)[x] *
probas_[nm];
403 if (!hasParameter(variable))
404 throw ParameterNotFoundException(
"DRHomogeneousTreeLikelihood::getFirstOrderDerivative().", variable);
405 if (getRateDistributionParameters().hasParameter(variable))
407 throw Exception(
"Derivatives respective to rate distribution parameters are not implemented.");
409 if (getSubstitutionModelParameters().hasParameter(variable))
411 throw Exception(
"Derivatives respective to substitution model parameters are not implemented.");
419 unsigned int brI = TextTools::to<unsigned int>(variable.substr(5));
420 const Node* branch = nodes_[brI];
421 vector< Vdouble*> _vdLikelihoods_branch;
422 for (
size_t i = 0; i < treeLikelihoodsContainer_.size(); i++)
424 _vdLikelihoods_branch.push_back(&treeLikelihoodsContainer_[i]->likelihoodData_->getDLikelihoodArray(branch->getId()));
429 const vector<unsigned int> * w = &likelihoodData_->getWeights();
430 for (
unsigned int i = 0; i < nbDistinctSites_; i++)
433 for (
unsigned int j = 0; j < treeLikelihoodsContainer_.size(); j++)
435 x += (*_vdLikelihoods_branch[j])[i] * probas_[j];
467 if (!hasParameter(variable))
468 throw ParameterNotFoundException(
"DRHomogeneousTreeLikelihood::getFirstOrderDerivative().", variable);
469 if (getRateDistributionParameters().hasParameter(variable))
471 throw Exception(
"Derivatives respective to rate distribution parameters are not implemented.");
473 if (getSubstitutionModelParameters().hasParameter(variable))
475 throw Exception(
"Derivatives respective to substitution model parameters are not implemented.");
483 unsigned int brI = TextTools::to<unsigned int>(variable.substr(5));
484 const Node* branch = nodes_[brI];
485 vector< Vdouble*> _vdLikelihoods_branch, _vd2Likelihoods_branch;
486 for (
unsigned int i = 0; i < treeLikelihoodsContainer_.size(); i++)
488 _vdLikelihoods_branch.push_back(&treeLikelihoodsContainer_[i]->likelihoodData_->getDLikelihoodArray(branch->getId()));
489 _vd2Likelihoods_branch.push_back(&treeLikelihoodsContainer_[i]->likelihoodData_->getD2LikelihoodArray(branch->getId()));
494 const vector<unsigned int> * w = &likelihoodData_->getWeights();
495 for (
unsigned int i = 0; i < nbDistinctSites_; i++)
499 for (
unsigned int j = 0; j < treeLikelihoodsContainer_.size(); j++)
501 x += (*_vdLikelihoods_branch[j])[i] * probas_[j];
503 for (
unsigned int j = 0; j < treeLikelihoodsContainer_.size(); j++)
505 x2 += (*_vd2Likelihoods_branch[j])[i] * probas_[j];
508 d += (*w)[i] * (x2 - pow(x, 2));
DRHomogeneousMixedTreeLikelihood & operator=(const DRHomogeneousMixedTreeLikelihood &lik)
DRHomogeneousTreeLikelihood & operator=(const DRHomogeneousTreeLikelihood &lik)
double getFirstOrderDerivative(const std::string &variable) const
Interface for all substitution models.
virtual void computeLikelihoodAtNode_(const Node *node, VVVdouble &likelihoodArray, const Node *sonNode=0) const
virtual double getNProbability(size_t i) const =0
Returns the probability of a specific model from the mixture.
double getLikelihoodForASite(size_t site) const
Get the likelihood for a site.
DRASDRTreeLikelihoodData * likelihoodData_
Vdouble & getRootRateSiteLikelihoodArray()
virtual void computeTreeD2LikelihoodAtNode(const Node *)
void computeTreeLikelihood()
virtual void computeSubtreeLikelihoodPostfix(const Node *node)
Compute the likelihood for a subtree defined by the Tree::Node node.
std::vector< DRHomogeneousTreeLikelihood * > treeLikelihoodsContainer_
const std::vector< unsigned int > & getWeights() const
Interface for phylogenetic tree objects.
SubstitutionModel * model_
double getLogLikelihoodForASite(size_t site) const
Get the logarithm of the likelihood for a site.
double getLogLikelihoodForASiteForARateClass(size_t site, size_t rateClass) const
Get the logarithm of the likelihood for a site knowing its rate class.
void initialize()
Init the likelihood object.
double getSecondOrderDerivative(const std::string &variable) const
void setData(const SiteContainer &sites)
Set the dataset for which the likelihood must be evaluated.
DRHomogeneousMixedTreeLikelihood(const Tree &tree, SubstitutionModel *model, DiscreteDistribution *rDist, bool checkRooted=true, bool verbose=true, bool rootArray=false)
Build a new DRHomogeneousMixedTreeLikelihood object without data.
void initialize()
Init the likelihood object.
virtual void computeTreeDLikelihoods()
virtual void computeSubtreeLikelihoodPrefix(const Node *node)
virtual void computeTreeD2Likelihoods()
double getLikelihoodForASiteForARateClass(size_t site, size_t rateClass) const
Get the likelihood for a site knowing its rate class.
This class implements the likelihood computation for a tree using the double-recursive algorithm...
double getLikelihood() const
Get the likelihood for the whole dataset.
virtual const std::vector< double > & getProbabilities() const =0
void fireParameterChanged(const ParameterList ¶ms)
virtual void displayLikelihood(const Node *node)
This method is mainly for debugging purpose.
void resetLikelihoodArrays(const Node *node)
The phylogenetic node class.
virtual void computeRootLikelihood()
virtual ~DRHomogeneousMixedTreeLikelihood()
void setData(const SiteContainer &sites)
Set the dataset for which the likelihood must be evaluated.
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.
std::vector< double > probas_
virtual size_t getNumberOfModels() const =0
virtual void applyParameters()
All parameters are stored in a parameter list. This function apply these parameters to the substituti...
A class to compute the average of several DRHomogeneousTreeLikelihood defined from a Mixed Substituti...
double getLogLikelihood() const
Get the logarithm of the likelihood for the whole dataset.
virtual const SubstitutionModel * getNModel(size_t i) const =0
Returns a specific model from the mixture.
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...
Interface for Substitution models, defined as a mixture of "simple" substitution models.