bpp-phyl  2.2.0
RHomogeneousTreeLikelihood.cpp
Go to the documentation of this file.
1 //
2 // File: RHomogeneousTreeLikelihood.cpp
3 // Created by: Julien Dutheil
4 // Created on: Fri Oct 17 18:14:51 2003
5 //
6 
7 /*
8  Copyright or © or Copr. Bio++ Development Team, (November 16, 2004)
9 
10  This software is a computer program whose purpose is to provide classes
11  for phylogenetic data analysis.
12 
13  This software is governed by the CeCILL license under French law and
14  abiding by the rules of distribution of free software. You can use,
15  modify and/ or redistribute the software under the terms of the CeCILL
16  license as circulated by CEA, CNRS and INRIA at the following URL
17  "http://www.cecill.info".
18 
19  As a counterpart to the access to the source code and rights to copy,
20  modify and redistribute granted by the license, users are provided only
21  with a limited warranty and the software's author, the holder of the
22  economic rights, and the successive licensors have only limited
23  liability.
24 
25  In this respect, the user's attention is drawn to the risks associated
26  with loading, using, modifying and/or developing or reproducing the
27  software by the user in light of its specific status of free software,
28  that may mean that it is complicated to manipulate, and that also
29  therefore means that it is reserved for developers and experienced
30  professionals having in-depth computer knowledge. Users are therefore
31  encouraged to load and test the software's suitability as regards their
32  requirements in conditions enabling the security of their systems and/or
33  data to be ensured and, more generally, to use and operate it in the
34  same conditions as regards security.
35 
36  The fact that you are presently reading this means that you have had
37  knowledge of the CeCILL license and that you accept its terms.
38  */
39 
41 #include "../PatternTools.h"
42 
43 #include <Bpp/Text/TextTools.h>
44 #include <Bpp/App/ApplicationTools.h>
45 
46 using namespace bpp;
47 
48 // From the STL:
49 #include <iostream>
50 
51 using namespace std;
52 
53 /******************************************************************************/
54 
56  const Tree& tree,
57  SubstitutionModel* model,
58  DiscreteDistribution* rDist,
59  bool checkRooted,
60  bool verbose,
61  bool usePatterns)
62 throw (Exception) :
63  AbstractHomogeneousTreeLikelihood(tree, model, rDist, checkRooted, verbose),
64  likelihoodData_(0),
65  minusLogLik_(-1.)
66 {
67  init_(usePatterns);
68 }
69 
70 /******************************************************************************/
71 
73  const Tree& tree,
74  const SiteContainer& data,
75  SubstitutionModel* model,
76  DiscreteDistribution* rDist,
77  bool checkRooted,
78  bool verbose,
79  bool usePatterns)
80 throw (Exception) :
81  AbstractHomogeneousTreeLikelihood(tree, model, rDist, checkRooted, verbose),
82  likelihoodData_(0),
83  minusLogLik_(-1.)
84 {
85  init_(usePatterns);
86  setData(data);
87 }
88 
89 /******************************************************************************/
90 
91 void RHomogeneousTreeLikelihood::init_(bool usePatterns) throw (Exception)
92 {
93  likelihoodData_ = new DRASRTreeLikelihoodData(
94  tree_,
95  rateDistribution_->getNumberOfCategories(),
96  usePatterns);
97 }
98 
99 /******************************************************************************/
100 
102  const RHomogeneousTreeLikelihood& lik) :
104  likelihoodData_(0),
105  minusLogLik_(lik.minusLogLik_)
106 {
109 }
110 
111 /******************************************************************************/
112 
114  const RHomogeneousTreeLikelihood& lik)
115 {
117  if (likelihoodData_) delete likelihoodData_;
121  return *this;
122 }
123 
124 /******************************************************************************/
125 
127 {
128  delete likelihoodData_;
129 }
130 
131 /******************************************************************************/
132 
133 void RHomogeneousTreeLikelihood::setData(const SiteContainer& sites) throw (Exception)
134 {
135  if (data_) delete data_;
136  data_ = PatternTools::getSequenceSubset(sites, *tree_->getRootNode());
137  if (verbose_) ApplicationTools::displayTask("Initializing data structure");
138  likelihoodData_->initLikelihoods(*data_, *model_);
139  if (verbose_) ApplicationTools::displayTaskDone();
140 
141  nbSites_ = likelihoodData_->getNumberOfSites();
142  nbDistinctSites_ = likelihoodData_->getNumberOfDistinctSites();
143  nbStates_ = likelihoodData_->getNumberOfStates();
144 
145  if (verbose_) ApplicationTools::displayResult("Number of distinct sites",
146  TextTools::toString(nbDistinctSites_));
147  initialized_ = false;
148 }
149 
150 /******************************************************************************/
151 
153 {
154  double l = 1.;
155  for (size_t i = 0; i < nbSites_; i++)
156  {
157  l *= getLikelihoodForASite(i);
158  }
159  return l;
160 }
161 
162 /******************************************************************************/
163 
165 {
166  double ll = 0;
167  vector<double> la(nbSites_);
168  for (size_t i = 0; i < nbSites_; i++)
169  {
170  la[i] = getLogLikelihoodForASite(i);
171  }
172  sort(la.begin(), la.end());
173  for (size_t i = nbSites_; i > 0; i--)
174  {
175  ll += la[i - 1];
176  }
177  return ll;
178 }
179 
180 /******************************************************************************/
181 
183 {
184  double l = 0;
185  for (size_t i = 0; i < nbClasses_; i++)
186  {
187  l += getLikelihoodForASiteForARateClass(site, i) * rateDistribution_->getProbability(i);
188  }
189  return l;
190 }
191 
192 /******************************************************************************/
193 
195 {
196  double l = 0;
197  for (size_t i = 0; i < nbClasses_; i++)
198  {
199  double li = getLikelihoodForASiteForARateClass(site, i) * rateDistribution_->getProbability(i);
200  if (li > 0) l+= li; //Corrects for numerical instabilities leading to slightly negative likelihoods
201  }
202  return log(l);
203 }
204 
205 /******************************************************************************/
206 
207 double RHomogeneousTreeLikelihood::getLikelihoodForASiteForARateClass(size_t site, size_t rateClass) const
208 {
209  double l = 0;
210  Vdouble* la = &likelihoodData_->getLikelihoodArray(tree_->getRootNode()->getId())[likelihoodData_->getRootArrayPosition(site)][rateClass];
211  for (size_t i = 0; i < nbStates_; i++)
212  {
213  //cout << (*la)[i] << "\t" << rootFreqs_[i] << endl;
214  double li = (*la)[i] * rootFreqs_[i];
215  if (li > 0) l+= li; //Corrects for numerical instabilities leading to slightly negative likelihoods
216  }
217  return l;
218 }
219 
220 /******************************************************************************/
221 
222 double RHomogeneousTreeLikelihood::getLogLikelihoodForASiteForARateClass(size_t site, size_t rateClass) const
223 {
224  double l = 0;
225  Vdouble* la = &likelihoodData_->getLikelihoodArray(tree_->getRootNode()->getId())[likelihoodData_->getRootArrayPosition(site)][rateClass];
226  for (size_t i = 0; i < nbStates_; i++)
227  {
228  l += (*la)[i] * rootFreqs_[i];
229  }
230  //if(l <= 0.) cerr << "WARNING!!! Negative likelihood." << endl;
231  return log(l);
232 }
233 
234 /******************************************************************************/
235 
236 double RHomogeneousTreeLikelihood::getLikelihoodForASiteForARateClassForAState(size_t site, size_t rateClass, int state) const
237 {
238  return likelihoodData_->getLikelihoodArray(tree_->getRootNode()->getId())[likelihoodData_->getRootArrayPosition(site)][rateClass][static_cast<size_t>(state)];
239 }
240 
241 /******************************************************************************/
242 
243 double RHomogeneousTreeLikelihood::getLogLikelihoodForASiteForARateClassForAState(size_t site, size_t rateClass, int state) const
244 {
245  return log(likelihoodData_->getLikelihoodArray(tree_->getRootNode()->getId())[likelihoodData_->getRootArrayPosition(site)][rateClass][static_cast<size_t>(state)]);
246 }
247 
248 /******************************************************************************/
249 
250 void RHomogeneousTreeLikelihood::setParameters(const ParameterList& parameters)
251 throw (ParameterNotFoundException, ConstraintException)
252 {
253  setParametersValues(parameters);
254 }
255 
256 /******************************************************************************/
257 
258 void RHomogeneousTreeLikelihood::fireParameterChanged(const ParameterList& params)
259 {
260  applyParameters();
261 
262  if (rateDistribution_->getParameters().getCommonParametersWith(params).size() > 0
263  || model_->getParameters().getCommonParametersWith(params).size() > 0)
264  {
265  //Rate parameter changed, need to recompute all probs:
267  }
268  else if (params.size() > 0)
269  {
270  //We may save some computations:
271  for (size_t i = 0; i < params.size(); i++)
272  {
273  string s = params[i].getName();
274  if (s.substr(0, 5) == "BrLen")
275  {
276  //Branch length parameter:
277  computeTransitionProbabilitiesForNode(nodes_[TextTools::to<size_t>(s.substr(5))]);
278  }
279  }
281  }
282 
284 
286 }
287 
288 /******************************************************************************/
289 
291 throw (Exception)
292 {
293  if (!isInitialized()) throw Exception("RHomogeneousTreeLikelihood::getValue(). Instance is not initialized.");
294  return minusLogLik_;
295 }
296 
297 /******************************************************************************
298 * First Order Derivatives *
299 ******************************************************************************/
300 
302  size_t site,
303  size_t rateClass) const
304 {
305  double dl = 0;
306  Vdouble* dla = &likelihoodData_->getDLikelihoodArray(tree_->getRootNode()->getId())[likelihoodData_->getRootArrayPosition(site)][rateClass];
307  for (size_t i = 0; i < nbStates_; i++)
308  {
309  dl += (*dla)[i] * rootFreqs_[i];
310  }
311  return dl;
312 }
313 
314 /******************************************************************************/
315 
317 {
318  // Derivative of the sum is the sum of derivatives:
319  double dl = 0;
320  for (size_t i = 0; i < nbClasses_; i++)
321  {
322  dl += getDLikelihoodForASiteForARateClass(site, i) * rateDistribution_->getProbability(i);
323  }
324  return dl;
325 }
326 
327 /******************************************************************************/
328 
330 {
331  // d(f(g(x)))/dx = dg(x)/dx . df(g(x))/dg :
332  return getDLikelihoodForASite(site) / getLikelihoodForASite(site);
333 }
334 
335 /******************************************************************************/
336 
338 {
339  // Derivative of the sum is the sum of derivatives:
340  double dl = 0;
341  for (size_t i = 0; i < nbSites_; i++)
342  {
343  dl += getDLogLikelihoodForASite(i);
344  }
345  return dl;
346 }
347 
348 /******************************************************************************/
349 
350 double RHomogeneousTreeLikelihood::getFirstOrderDerivative(const string& variable) const
351 throw (Exception)
352 {
353  if (!hasParameter(variable))
354  throw ParameterNotFoundException("RHomogeneousTreeLikelihood::getFirstOrderDerivative().", variable);
355  if (getRateDistributionParameters().hasParameter(variable))
356  {
357  throw Exception("Derivatives respective to rate distribution parameter are not implemented.");
358  }
359  if (getSubstitutionModelParameters().hasParameter(variable))
360  {
361  throw Exception("Derivatives respective to substitution model parameters are not implemented.");
362  }
363 
364  const_cast<RHomogeneousTreeLikelihood*>(this)->computeTreeDLikelihood(variable);
365  return -getDLogLikelihood();
366 }
367 
368 /******************************************************************************/
369 
371 {
372  // Get the node with the branch whose length must be derivated:
373  size_t brI = TextTools::to<size_t>(variable.substr(5));
374  const Node* branch = nodes_[brI];
375  const Node* father = branch->getFather();
376  VVVdouble* _dLikelihoods_father = &likelihoodData_->getDLikelihoodArray(father->getId());
377 
378  // Compute dLikelihoods array for the father node.
379  // Fist initialize to 1:
380  size_t nbSites = _dLikelihoods_father->size();
381  for (size_t i = 0; i < nbSites; i++)
382  {
383  VVdouble* _dLikelihoods_father_i = &(*_dLikelihoods_father)[i];
384  for (size_t c = 0; c < nbClasses_; c++)
385  {
386  Vdouble* _dLikelihoods_father_i_c = &(*_dLikelihoods_father_i)[c];
387  for (size_t s = 0; s < nbStates_; s++)
388  {
389  (*_dLikelihoods_father_i_c)[s] = 1.;
390  }
391  }
392  }
393 
394  size_t nbNodes = father->getNumberOfSons();
395  for (size_t l = 0; l < nbNodes; l++)
396  {
397  const Node* son = father->getSon(l);
398 
399  vector<size_t> * _patternLinks_father_son = &likelihoodData_->getArrayPositions(father->getId(), son->getId());
400  VVVdouble* _likelihoods_son = &likelihoodData_->getLikelihoodArray(son->getId());
401 
402  if (son == branch)
403  {
404  VVVdouble* dpxy__son = &dpxy_[son->getId()];
405  for (size_t i = 0; i < nbSites; i++)
406  {
407  VVdouble* _likelihoods_son_i = &(*_likelihoods_son)[(*_patternLinks_father_son)[i]];
408  VVdouble* _dLikelihoods_father_i = &(*_dLikelihoods_father)[i];
409  for (size_t c = 0; c < nbClasses_; c++)
410  {
411  Vdouble* _likelihoods_son_i_c = &(*_likelihoods_son_i)[c];
412  Vdouble* _dLikelihoods_father_i_c = &(*_dLikelihoods_father_i)[c];
413  VVdouble* dpxy__son_c = &(*dpxy__son)[c];
414  for (size_t x = 0; x < nbStates_; x++)
415  {
416  double dl = 0;
417  Vdouble* dpxy__son_c_x = &(*dpxy__son_c)[x];
418  for (size_t y = 0; y < nbStates_; y++)
419  {
420  dl += (*dpxy__son_c_x)[y] * (*_likelihoods_son_i_c)[y];
421  }
422  (*_dLikelihoods_father_i_c)[x] *= dl;
423  }
424  }
425  }
426  }
427  else
428  {
429  VVVdouble* pxy__son = &pxy_[son->getId()];
430  for (size_t i = 0; i < nbSites; i++)
431  {
432  VVdouble* _likelihoods_son_i = &(*_likelihoods_son)[(*_patternLinks_father_son)[i]];
433  VVdouble* _dLikelihoods_father_i = &(*_dLikelihoods_father)[i];
434  for (size_t c = 0; c < nbClasses_; c++)
435  {
436  Vdouble* _likelihoods_son_i_c = &(*_likelihoods_son_i)[c];
437  Vdouble* _dLikelihoods_father_i_c = &(*_dLikelihoods_father_i)[c];
438  VVdouble* pxy__son_c = &(*pxy__son)[c];
439  for (size_t x = 0; x < nbStates_; x++)
440  {
441  double dl = 0;
442  Vdouble* pxy__son_c_x = &(*pxy__son_c)[x];
443  for (size_t y = 0; y < nbStates_; y++)
444  {
445  dl += (*pxy__son_c_x)[y] * (*_likelihoods_son_i_c)[y];
446  }
447  (*_dLikelihoods_father_i_c)[x] *= dl;
448  }
449  }
450  }
451  }
452  }
453 
454  // Now we go down the tree toward the root node:
456 }
457 
458 /******************************************************************************/
459 
461 {
462  const Node* father = node->getFather();
463  // We assume that the _dLikelihoods array has been filled for the current node 'node'.
464  // We will evaluate the array for the father node.
465  if (father == NULL) return; // We reached the root!
466 
467  // Compute dLikelihoods array for the father node.
468  // Fist initialize to 1:
469  VVVdouble* _dLikelihoods_father = &likelihoodData_->getDLikelihoodArray(father->getId());
470  size_t nbSites = _dLikelihoods_father->size();
471  for (size_t i = 0; i < nbSites; i++)
472  {
473  VVdouble* _dLikelihoods_father_i = &(*_dLikelihoods_father)[i];
474  for (size_t c = 0; c < nbClasses_; c++)
475  {
476  Vdouble* _dLikelihoods_father_i_c = &(*_dLikelihoods_father_i)[c];
477  for (size_t s = 0; s < nbStates_; s++)
478  {
479  (*_dLikelihoods_father_i_c)[s] = 1.;
480  }
481  }
482  }
483 
484  size_t nbNodes = father->getNumberOfSons();
485  for (size_t l = 0; l < nbNodes; l++)
486  {
487  const Node* son = father->getSon(l);
488 
489  VVVdouble* pxy__son = &pxy_[son->getId()];
490  vector<size_t> * _patternLinks_father_son = &likelihoodData_->getArrayPositions(father->getId(), son->getId());
491 
492  if (son == node)
493  {
494  VVVdouble* _dLikelihoods_son = &likelihoodData_->getDLikelihoodArray(son->getId());
495  for (size_t i = 0; i < nbSites; i++)
496  {
497  VVdouble* _dLikelihoods_son_i = &(*_dLikelihoods_son)[(*_patternLinks_father_son)[i]];
498  VVdouble* _dLikelihoods_father_i = &(*_dLikelihoods_father)[i];
499  for (size_t c = 0; c < nbClasses_; c++)
500  {
501  Vdouble* _dLikelihoods_son_i_c = &(*_dLikelihoods_son_i)[c];
502  Vdouble* _dLikelihoods_father_i_c = &(*_dLikelihoods_father_i)[c];
503  VVdouble* pxy__son_c = &(*pxy__son)[c];
504  for (size_t x = 0; x < nbStates_; x++)
505  {
506  double dl = 0;
507  Vdouble* pxy__son_c_x = &(*pxy__son_c)[x];
508  for (size_t y = 0; y < nbStates_; y++)
509  {
510  dl += (*pxy__son_c_x)[y] * (*_dLikelihoods_son_i_c)[y];
511  }
512  (*_dLikelihoods_father_i_c)[x] *= dl;
513  }
514  }
515  }
516  }
517  else
518  {
519  VVVdouble* _likelihoods_son = &likelihoodData_->getLikelihoodArray(son->getId());
520  for (size_t i = 0; i < nbSites; i++)
521  {
522  VVdouble* _likelihoods_son_i = &(*_likelihoods_son)[(*_patternLinks_father_son)[i]];
523  VVdouble* _dLikelihoods_father_i = &(*_dLikelihoods_father)[i];
524  for (size_t c = 0; c < nbClasses_; c++)
525  {
526  Vdouble* _likelihoods_son_i_c = &(*_likelihoods_son_i)[c];
527  Vdouble* _dLikelihoods_father_i_c = &(*_dLikelihoods_father_i)[c];
528  VVdouble* pxy__son_c = &(*pxy__son)[c];
529  for (size_t x = 0; x < nbStates_; x++)
530  {
531  double dl = 0;
532  Vdouble* pxy__son_c_x = &(*pxy__son_c)[x];
533  for (size_t y = 0; y < nbStates_; y++)
534  {
535  dl += (*pxy__son_c_x)[y] * (*_likelihoods_son_i_c)[y];
536  }
537  (*_dLikelihoods_father_i_c)[x] *= dl;
538  }
539  }
540  }
541  }
542  }
543 
544  //Next step: move toward grand father...
546 }
547 
548 /******************************************************************************
549 * Second Order Derivatives *
550 ******************************************************************************/
551 
553  size_t site,
554  size_t rateClass) const
555 {
556  double d2l = 0;
557  Vdouble* d2la = &likelihoodData_->getD2LikelihoodArray(tree_->getRootNode()->getId())[likelihoodData_->getRootArrayPosition(site)][rateClass];
558  for (size_t i = 0; i < nbStates_; i++)
559  {
560  d2l += (*d2la)[i] * rootFreqs_[i];
561  }
562  return d2l;
563 }
564 
565 /******************************************************************************/
566 
568 {
569  // Derivative of the sum is the sum of derivatives:
570  double d2l = 0;
571  for (size_t i = 0; i < nbClasses_; i++)
572  {
573  d2l += getD2LikelihoodForASiteForARateClass(site, i) * rateDistribution_->getProbability(i);
574  }
575  return d2l;
576 }
577 
578 /******************************************************************************/
579 
581 {
583  - pow( getDLikelihoodForASite(site) / getLikelihoodForASite(site), 2);
584 }
585 
586 /******************************************************************************/
587 
589 {
590  // Derivative of the sum is the sum of derivatives:
591  double dl = 0;
592  for (size_t i = 0; i < nbSites_; i++)
593  {
595  }
596  return dl;
597 }
598 
599 /******************************************************************************/
600 
601 double RHomogeneousTreeLikelihood::getSecondOrderDerivative(const string& variable) const
602 throw (Exception)
603 {
604  if (!hasParameter(variable))
605  throw ParameterNotFoundException("RHomogeneousTreeLikelihood::getSecondOrderDerivative().", variable);
606  if (getRateDistributionParameters().hasParameter(variable))
607  {
608  throw Exception("Derivatives respective to rate distribution parameter are not implemented.");
609  }
610  if (getSubstitutionModelParameters().hasParameter(variable))
611  {
612  throw Exception("Derivatives respective to substitution model parameters are not implemented.");
613  }
614 
615  const_cast<RHomogeneousTreeLikelihood*>(this)->computeTreeD2Likelihood(variable);
616  return -getD2LogLikelihood();
617 }
618 
619 /******************************************************************************/
620 
622 {
623  // Get the node with the branch whose length must be derivated:
624  size_t brI = TextTools::to<size_t>(variable.substr(5));
625  const Node* branch = nodes_[brI];
626  const Node* father = branch->getFather();
627 
628  // Compute dLikelihoods array for the father node.
629  // Fist initialize to 1:
630  VVVdouble* _d2Likelihoods_father = &likelihoodData_->getD2LikelihoodArray(father->getId());
631  size_t nbSites = _d2Likelihoods_father->size();
632  for (size_t i = 0; i < nbSites; i++)
633  {
634  VVdouble* _d2Likelihoods_father_i = &(*_d2Likelihoods_father)[i];
635  for (size_t c = 0; c < nbClasses_; c++)
636  {
637  Vdouble* _d2Likelihoods_father_i_c = &(*_d2Likelihoods_father_i)[c];
638  for (size_t s = 0; s < nbStates_; s++)
639  {
640  (*_d2Likelihoods_father_i_c)[s] = 1.;
641  }
642  }
643  }
644 
645  size_t nbNodes = father->getNumberOfSons();
646  for (size_t l = 0; l < nbNodes; l++)
647  {
648  const Node* son = father->getSon(l);
649 
650  vector<size_t> * _patternLinks_father_son = &likelihoodData_->getArrayPositions(father->getId(), son->getId());
651  VVVdouble* _likelihoods_son = &likelihoodData_->getLikelihoodArray(son->getId());
652 
653  if (son == branch)
654  {
655  VVVdouble* d2pxy__son = &d2pxy_[son->getId()];
656  for (size_t i = 0; i < nbSites; i++)
657  {
658  VVdouble* _likelihoods_son_i = &(*_likelihoods_son)[(*_patternLinks_father_son)[i]];
659  VVdouble* _d2Likelihoods_father_i = &(*_d2Likelihoods_father)[i];
660  for (size_t c = 0; c < nbClasses_; c++)
661  {
662  Vdouble* _likelihoods_son_i_c = &(*_likelihoods_son_i)[c];
663  Vdouble* _d2Likelihoods_father_i_c = &(*_d2Likelihoods_father_i)[c];
664  VVdouble* d2pxy__son_c = &(*d2pxy__son)[c];
665  for (size_t x = 0; x < nbStates_; x++)
666  {
667  double d2l = 0;
668  Vdouble* d2pxy__son_c_x = &(*d2pxy__son_c)[x];
669  for (size_t y = 0; y < nbStates_; y++)
670  {
671  d2l += (*d2pxy__son_c_x)[y] * (*_likelihoods_son_i_c)[y];
672  }
673  (*_d2Likelihoods_father_i_c)[x] *= d2l;
674  }
675  }
676  }
677  }
678  else
679  {
680  VVVdouble* pxy__son = &pxy_[son->getId()];
681  for (size_t i = 0; i < nbSites; i++)
682  {
683  VVdouble* _likelihoods_son_i = &(*_likelihoods_son)[(*_patternLinks_father_son)[i]];
684  VVdouble* _d2Likelihoods_father_i = &(*_d2Likelihoods_father)[i];
685  for (size_t c = 0; c < nbClasses_; c++)
686  {
687  Vdouble* _likelihoods_son_i_c = &(*_likelihoods_son_i)[c];
688  Vdouble* _d2Likelihoods_father_i_c = &(*_d2Likelihoods_father_i)[c];
689  VVdouble* pxy__son_c = &(*pxy__son)[c];
690  for (size_t x = 0; x < nbStates_; x++)
691  {
692  double d2l = 0;
693  Vdouble* pxy__son_c_x = &(*pxy__son_c)[x];
694  for (size_t y = 0; y < nbStates_; y++)
695  {
696  d2l += (*pxy__son_c_x)[y] * (*_likelihoods_son_i_c)[y];
697  }
698  (*_d2Likelihoods_father_i_c)[x] *= d2l;
699  }
700  }
701  }
702  }
703  }
704 
705  // Now we go down the tree toward the root node:
707 }
708 
709 /******************************************************************************/
710 
712 {
713  const Node* father = node->getFather();
714  // We assume that the _dLikelihoods array has been filled for the current node 'node'.
715  // We will evaluate the array for the father node.
716  if (father == NULL) return; // We reached the root!
717 
718  // Compute dLikelihoods array for the father node.
719  // Fist initialize to 1:
720  VVVdouble* _d2Likelihoods_father = &likelihoodData_->getD2LikelihoodArray(father->getId());
721  size_t nbSites = _d2Likelihoods_father->size();
722  for (size_t i = 0; i < nbSites; i++)
723  {
724  VVdouble* _d2Likelihoods_father_i = &(*_d2Likelihoods_father)[i];
725  for (size_t c = 0; c < nbClasses_; c++)
726  {
727  Vdouble* _d2Likelihoods_father_i_c = &(*_d2Likelihoods_father_i)[c];
728  for (size_t s = 0; s < nbStates_; s++)
729  {
730  (*_d2Likelihoods_father_i_c)[s] = 1.;
731  }
732  }
733  }
734 
735  size_t nbNodes = father->getNumberOfSons();
736  for (size_t l = 0; l < nbNodes; l++)
737  {
738  const Node* son = father->getSon(l);
739 
740  VVVdouble* pxy__son = &pxy_[son->getId()];
741  vector<size_t> * _patternLinks_father_son = &likelihoodData_->getArrayPositions(father->getId(), son->getId());
742 
743  if (son == node)
744  {
745  VVVdouble* _d2Likelihoods_son = &likelihoodData_->getD2LikelihoodArray(son->getId());
746  for (size_t i = 0; i < nbSites; i++)
747  {
748  VVdouble* _d2Likelihoods_son_i = &(*_d2Likelihoods_son)[(*_patternLinks_father_son)[i]];
749  VVdouble* _d2Likelihoods_father_i = &(*_d2Likelihoods_father)[i];
750  for (size_t c = 0; c < nbClasses_; c++)
751  {
752  Vdouble* _d2Likelihoods_son_i_c = &(*_d2Likelihoods_son_i)[c];
753  Vdouble* _d2Likelihoods_father_i_c = &(*_d2Likelihoods_father_i)[c];
754  VVdouble* pxy__son_c = &(*pxy__son)[c];
755  for (size_t x = 0; x < nbStates_; x++)
756  {
757  double d2l = 0;
758  Vdouble* pxy__son_c_x = &(*pxy__son_c)[x];
759  for (size_t y = 0; y < nbStates_; y++)
760  {
761  d2l += (*pxy__son_c_x)[y] * (*_d2Likelihoods_son_i_c)[y];
762  }
763  (*_d2Likelihoods_father_i_c)[x] *= d2l;
764  }
765  }
766  }
767  }
768  else
769  {
770  VVVdouble* _likelihoods_son = &likelihoodData_->getLikelihoodArray(son->getId());
771  for (size_t i = 0; i < nbSites; i++)
772  {
773  VVdouble* _likelihoods_son_i = &(*_likelihoods_son)[(*_patternLinks_father_son)[i]];
774  VVdouble* _d2Likelihoods_father_i = &(*_d2Likelihoods_father)[i];
775  for (size_t c = 0; c < nbClasses_; c++)
776  {
777  Vdouble* _likelihoods_son_i_c = &(*_likelihoods_son_i)[c];
778  Vdouble* _d2Likelihoods_father_i_c = &(*_d2Likelihoods_father_i)[c];
779  VVdouble* pxy__son_c = &(*pxy__son)[c];
780  for (size_t x = 0; x < nbStates_; x++)
781  {
782  double dl = 0;
783  Vdouble* pxy__son_c_x = &(*pxy__son_c)[x];
784  for (size_t y = 0; y < nbStates_; y++)
785  {
786  dl += (*pxy__son_c_x)[y] * (*_likelihoods_son_i_c)[y];
787  }
788  (*_d2Likelihoods_father_i_c)[x] *= dl;
789  }
790  }
791  }
792  }
793  }
794 
795  //Next step: move toward grand father...
797 }
798 
799 /******************************************************************************/
800 
802 {
803  computeSubtreeLikelihood(tree_->getRootNode());
804 }
805 
806 /******************************************************************************/
807 
809 {
810  if (node->isLeaf()) return;
811 
812  size_t nbSites = likelihoodData_->getLikelihoodArray(node->getId()).size();
813  size_t nbNodes = node->getNumberOfSons();
814 
815  // Must reset the likelihood array first (i.e. set all of them to 1):
816  VVVdouble* _likelihoods_node = &likelihoodData_->getLikelihoodArray(node->getId());
817  for (size_t i = 0; i < nbSites; i++)
818  {
819  //For each site in the sequence,
820  VVdouble* _likelihoods_node_i = &(*_likelihoods_node)[i];
821  for (size_t c = 0; c < nbClasses_; c++)
822  {
823  //For each rate classe,
824  Vdouble* _likelihoods_node_i_c = &(*_likelihoods_node_i)[c];
825  for (size_t x = 0; x < nbStates_; x++)
826  {
827  //For each initial state,
828  (*_likelihoods_node_i_c)[x] = 1.;
829  }
830  }
831  }
832 
833  for (size_t l = 0; l < nbNodes; l++)
834  {
835  //For each son node,
836 
837  const Node* son = node->getSon(l);
838 
839  computeSubtreeLikelihood(son); //Recursive method:
840 
841  VVVdouble* pxy__son = &pxy_[son->getId()];
842  vector<size_t> * _patternLinks_node_son = &likelihoodData_->getArrayPositions(node->getId(), son->getId());
843  VVVdouble* _likelihoods_son = &likelihoodData_->getLikelihoodArray(son->getId());
844 
845  for (size_t i = 0; i < nbSites; i++)
846  {
847  //For each site in the sequence,
848  VVdouble* _likelihoods_son_i = &(*_likelihoods_son)[(*_patternLinks_node_son)[i]];
849  VVdouble* _likelihoods_node_i = &(*_likelihoods_node)[i];
850  for (size_t c = 0; c < nbClasses_; c++)
851  {
852  //For each rate classe,
853  Vdouble* _likelihoods_son_i_c = &(*_likelihoods_son_i)[c];
854  Vdouble* _likelihoods_node_i_c = &(*_likelihoods_node_i)[c];
855  VVdouble* pxy__son_c = &(*pxy__son)[c];
856  for (size_t x = 0; x < nbStates_; x++)
857  {
858  //For each initial state,
859  Vdouble* pxy__son_c_x = &(*pxy__son_c)[x];
860  double likelihood = 0;
861  for (size_t y = 0; y < nbStates_; y++)
862  {
863  likelihood += (*pxy__son_c_x)[y] * (*_likelihoods_son_i_c)[y];
864  }
865  (*_likelihoods_node_i_c)[x] *= likelihood;
866  }
867  }
868  }
869  }
870 }
871 
872 /******************************************************************************/
873 
875 {
876  cout << "Likelihoods at node " << node->getName() << ": " << endl;
878  cout << " ***" << endl;
879 }
880 
881 /*******************************************************************************/
882 
virtual void computeDownSubtreeDLikelihood(const Node *)
Interface for all substitution models.
virtual void computeAllTransitionProbabilities()
Fill the pxy_, dpxy_ and d2pxy_ arrays for all nodes.
RHomogeneousTreeLikelihood(const Tree &tree, SubstitutionModel *model, DiscreteDistribution *rDist, bool checkRooted=true, bool verbose=true, bool usePatterns=true)
Build a new RHomogeneousTreeLikelihood object without data.
double getLogLikelihoodForASite(size_t site) const
Get the logarithm of the likelihood for a site.
virtual std::string getName() const
Get the name associated to this node, if there is one, otherwise throw a NodeException.
Definition: Node.h:236
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...
size_t getRootArrayPosition(size_t currentPosition) const
DRASRTreeLikelihoodData * likelihoodData_
AbstractHomogeneousTreeLikelihood & operator=(const AbstractHomogeneousTreeLikelihood &lik)
Assignation operator.
void setTree(const TreeTemplate< Node > *tree)
Set the tree associated to the data.
virtual double getD2LogLikelihoodForASite(size_t site) const
virtual void computeDownSubtreeD2Likelihood(const Node *)
virtual double getD2LikelihoodForASiteForARateClass(size_t site, size_t rateClass) const
virtual void computeTreeD2Likelihood(const std::string &variable)
virtual const Node * getSon(size_t pos) const
Definition: Node.h:395
STL namespace.
discrete Rate Across Sites, (simple) Recursive likelihood data structure.
double getLikelihood() const
Get the likelihood for the whole dataset.
virtual const Node * getFather() const
Get the father of this node is there is one.
Definition: Node.h:339
Interface for phylogenetic tree objects.
Definition: Tree.h:148
static void displayLikelihoodArray(const VVVdouble &likelihoodArray)
Print the likelihood array to terminal (debugging tool).
This class implement the &#39;traditional&#39; way of computing likelihood for a tree.
virtual const Vdouble & getFrequencies() const =0
virtual bool isLeaf() const
Definition: Node.h:692
virtual double getDLogLikelihoodForASite(size_t site) const
RHomogeneousTreeLikelihood & operator=(const RHomogeneousTreeLikelihood &lik)
VVVdouble & getDLikelihoodArray(int nodeId)
void fireParameterChanged(const ParameterList &params)
virtual int getId() const
Get the node&#39;s id.
Definition: Node.h:203
double getLogLikelihood() const
Get the logarithm of the likelihood for the whole dataset.
static SiteContainer * getSequenceSubset(const SiteContainer &sequenceSet, const Node &node)
Extract the sequences corresponding to a given subtree.
virtual double getD2LikelihoodForASite(size_t site) const
The phylogenetic node class.
Definition: Node.h:90
double getLikelihoodForASite(size_t site) const
Get the likelihood for a site.
void setParameters(const ParameterList &parameters)
Implements the Function interface.
virtual void computeTransitionProbabilitiesForNode(const Node *node)
Fill the pxy_, dpxy_ and d2pxy_ arrays for one node.
double getLogLikelihoodForASiteForARateClass(size_t site, size_t rateClass) const
Get the logarithm of the likelihood for a site knowing its rate class.
DRASRTreeLikelihoodData * clone() const
virtual void computeTreeDLikelihood(const std::string &variable)
std::vector< Node * > nodes_
Pointer toward all nodes in the tree.
virtual void computeSubtreeLikelihood(const Node *node)
Compute the likelihood for a subtree defined by the Tree::Node node.
virtual double getDLikelihoodForASiteForARateClass(size_t site, size_t rateClass) const
double getSecondOrderDerivative(const std::string &variable) const
virtual size_t getNumberOfSons() const
Definition: Node.h:388
virtual double getDLikelihoodForASite(size_t site) const
double getFirstOrderDerivative(const std::string &variable) const
virtual void displayLikelihood(const Node *node)
This method is mainly for debugging purpose.
const std::vector< size_t > & getArrayPositions(int parentId, int sonId) const
virtual void applyParameters()
All parameters are stored in a parameter list. This function apply these parameters to the substituti...
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.
VVVdouble & getLikelihoodArray(int nodeId)
Partial implementation for homogeneous model of the TreeLikelihood interface.
VVVdouble & getD2LikelihoodArray(int nodeId)
void init_(bool usePatterns)
Method called by constructors.
double getLikelihoodForASiteForARateClass(size_t site, size_t rateClass) const
Get the likelihood for a site knowing its rate class.
void setData(const SiteContainer &sites)
Set the dataset for which the likelihood must be evaluated.