bpp-phyl3 3.0.0
RHomogeneousTreeLikelihood.cpp
Go to the documentation of this file.
1// SPDX-FileCopyrightText: The Bio++ Development Group
2//
3// SPDX-License-Identifier: CECILL-2.1
4
7
8#include "../../PatternTools.h"
10
11using namespace bpp;
12
13// From the STL:
14#include <iostream>
15
16using namespace std;
17
18/******************************************************************************/
19
20RHomogeneousTreeLikelihood::RHomogeneousTreeLikelihood(
21 const Tree& tree,
22 std::shared_ptr<TransitionModelInterface> model,
23 std::shared_ptr<DiscreteDistributionInterface> rDist,
24 bool checkRooted,
25 bool verbose,
26 bool usePatterns) :
27 AbstractHomogeneousTreeLikelihood(tree, model, rDist, checkRooted, verbose),
28 likelihoodData_(),
29 minusLogLik_(-1.)
30{
31 init_(usePatterns);
32}
33
34/******************************************************************************/
35
37 const Tree& tree,
38 const AlignmentDataInterface& data,
39 std::shared_ptr<TransitionModelInterface> model,
40 std::shared_ptr<DiscreteDistributionInterface> rDist,
41 bool checkRooted,
42 bool verbose,
43 bool usePatterns) :
44 AbstractHomogeneousTreeLikelihood(tree, model, rDist, checkRooted, verbose),
45 likelihoodData_(),
46 minusLogLik_(-1.)
47{
48 init_(usePatterns);
50}
51
52/******************************************************************************/
53
55{
57 tree_,
58 rateDistribution_->getNumberOfCategories(),
59 usePatterns);
60}
61
62/******************************************************************************/
63
65 const RHomogeneousTreeLikelihood& lik) :
67 likelihoodData_(0),
68 minusLogLik_(lik.minusLogLik_)
69{
72}
73
74/******************************************************************************/
75
78{
81 delete likelihoodData_;
85 return *this;
86}
87
88/******************************************************************************/
89
91{
92 delete likelihoodData_;
93}
94
95/******************************************************************************/
96
98{
99 data_ = PatternTools::getSequenceSubset(sites, *tree_->getRootNode());
100 if (verbose_)
101 ApplicationTools::displayTask("Initializing data structure");
103 if (verbose_)
105
109
110 if (verbose_)
111 ApplicationTools::displayResult("Number of distinct sites",
113 initialized_ = false;
114}
115
116/******************************************************************************/
117
119{
120 double l = 1.;
121 for (size_t i = 0; i < nbSites_; i++)
122 {
123 l *= getLikelihoodForASite(i);
124 }
125 return l;
126}
127
128/******************************************************************************/
129
131{
132 double ll = 0;
133 vector<double> la(nbSites_);
134 for (size_t i = 0; i < nbSites_; i++)
135 {
136 la[i] = getLogLikelihoodForASite(i);
137 }
138 sort(la.begin(), la.end());
139 for (size_t i = nbSites_; i > 0; i--)
140 {
141 ll += la[i - 1];
142 }
143 return ll;
144}
145
146/******************************************************************************/
147
149{
150 double l = 0;
151 for (size_t i = 0; i < nbClasses_; i++)
152 {
153 l += getLikelihoodForASiteForARateClass(site, i) * rateDistribution_->getProbability(i);
154 }
155 return l;
156}
157
158/******************************************************************************/
159
161{
162 double l = 0;
163 for (size_t i = 0; i < nbClasses_; i++)
164 {
165 double li = getLikelihoodForASiteForARateClass(site, i) * rateDistribution_->getProbability(i);
166 if (li > 0)
167 l += li; // Corrects for numerical instabilities leading to slightly negative likelihoods
168 }
169 return log(l);
170}
171
172/******************************************************************************/
173
174double RHomogeneousTreeLikelihood::getLikelihoodForASiteForARateClass(size_t site, size_t rateClass) const
175{
176 double l = 0;
177 Vdouble* la = &likelihoodData_->getLikelihoodArray(tree_->getRootNode()->getId())[likelihoodData_->getRootArrayPosition(site)][rateClass];
178 for (size_t i = 0; i < nbStates_; i++)
179 {
180 // cout << (*la)[i] << "\t" << rootFreqs_[i] << endl;
181 double li = (*la)[i] * rootFreqs_[i];
182 if (li > 0)
183 l += li; // Corrects for numerical instabilities leading to slightly negative likelihoods
184 }
185 return l;
186}
187
188/******************************************************************************/
189
191{
192 double l = 0;
193 Vdouble* la = &likelihoodData_->getLikelihoodArray(tree_->getRootNode()->getId())[likelihoodData_->getRootArrayPosition(site)][rateClass];
194 for (size_t i = 0; i < nbStates_; i++)
195 {
196 l += (*la)[i] * rootFreqs_[i];
197 }
198 // if(l <= 0.) cerr << "WARNING!!! Negative likelihood." << endl;
199 return log(l);
200}
201
202/******************************************************************************/
203
204double RHomogeneousTreeLikelihood::getLikelihoodForASiteForARateClassForAState(size_t site, size_t rateClass, int state) const
205{
206 return likelihoodData_->getLikelihoodArray(tree_->getRootNode()->getId())[likelihoodData_->getRootArrayPosition(site)][rateClass][static_cast<size_t>(state)];
207}
208
209/******************************************************************************/
210
211double RHomogeneousTreeLikelihood::getLogLikelihoodForASiteForARateClassForAState(size_t site, size_t rateClass, int state) const
212{
213 return log(likelihoodData_->getLikelihoodArray(tree_->getRootNode()->getId())[likelihoodData_->getRootArrayPosition(site)][rateClass][static_cast<size_t>(state)]);
214}
215
216/******************************************************************************/
217
219{
220 setParametersValues(parameters);
221}
222
223/******************************************************************************/
224
226{
228
229 if (rateDistribution_->getParameters().getCommonParametersWith(params).size() > 0
230 || model_->getParameters().getCommonParametersWith(params).size() > 0)
231 {
232 // Rate parameter changed, need to recompute all probs:
234 }
235 else if (params.size() > 0)
236 {
237 // We may save some computations:
238 for (size_t i = 0; i < params.size(); i++)
239 {
240 string s = params[i].getName();
241 if (s.substr(0, 5) == "BrLen")
242 {
243 // Branch length parameter:
244 computeTransitionProbabilitiesForNode(nodes_[TextTools::to<size_t>(s.substr(5))]);
245 }
246 }
247 rootFreqs_ = model_->getFrequencies();
248 }
249
251
253}
254
255/******************************************************************************/
256
258{
259 if (!isInitialized())
260 throw Exception("RHomogeneousTreeLikelihood::getValue(). Instance is not initialized.");
261 return minusLogLik_;
262}
263
264/******************************************************************************
265 * First Order Derivatives *
266 ******************************************************************************/
268 size_t site,
269 size_t rateClass) const
270{
271 double dl = 0;
272 Vdouble* dla = &likelihoodData_->getDLikelihoodArray(tree_->getRootNode()->getId())[likelihoodData_->getRootArrayPosition(site)][rateClass];
273 for (size_t i = 0; i < nbStates_; i++)
274 {
275 dl += (*dla)[i] * rootFreqs_[i];
276 }
277 return dl;
278}
279
280/******************************************************************************/
281
283{
284 // Derivative of the sum is the sum of derivatives:
285 double dl = 0;
286 for (size_t i = 0; i < nbClasses_; i++)
287 {
288 dl += getDLikelihoodForASiteForARateClass(site, i) * rateDistribution_->getProbability(i);
289 }
290 return dl;
291}
292
293/******************************************************************************/
294
296{
297 // d(f(g(x)))/dx = dg(x)/dx . df(g(x))/dg :
299}
300
301/******************************************************************************/
302
304{
305 // Derivative of the sum is the sum of derivatives:
306 double dl = 0;
307 for (size_t i = 0; i < nbSites_; i++)
308 {
310 }
311
312 return dl;
313}
314
315/******************************************************************************/
316
317double RHomogeneousTreeLikelihood::getFirstOrderDerivative(const string& variable) const
318{
319 if (!hasParameter(variable))
320 throw ParameterNotFoundException("RHomogeneousTreeLikelihood::getFirstOrderDerivative().", variable);
322 {
323 throw Exception("Derivatives respective to rate distribution parameter are not implemented.");
324 }
326 {
327 throw Exception("Derivatives respective to substitution model parameters are not implemented.");
328 }
329
330 const_cast<RHomogeneousTreeLikelihood*>(this)->computeTreeDLikelihood(variable);
331
332 return -getDLogLikelihood();
333}
334
335/******************************************************************************/
336
338{
339 // Get the node with the branch whose length must be derivated:
340 size_t brI = TextTools::to<size_t>(variable.substr(5));
341 const Node* branch = nodes_[brI];
342 const Node* father = branch->getFather();
343 VVVdouble* _dLikelihoods_father = &likelihoodData_->getDLikelihoodArray(father->getId());
344
345 // Compute dLikelihoods array for the father node.
346 // Fist initialize to 1:
347 size_t nbSites = _dLikelihoods_father->size();
348 for (size_t i = 0; i < nbSites; i++)
349 {
350 VVdouble* _dLikelihoods_father_i = &(*_dLikelihoods_father)[i];
351 for (size_t c = 0; c < nbClasses_; c++)
352 {
353 Vdouble* _dLikelihoods_father_i_c = &(*_dLikelihoods_father_i)[c];
354 for (size_t s = 0; s < nbStates_; s++)
355 {
356 (*_dLikelihoods_father_i_c)[s] = 1.;
357 }
358 }
359 }
360
361 size_t nbNodes = father->getNumberOfSons();
362 for (size_t l = 0; l < nbNodes; l++)
363 {
364 const Node* son = father->getSon(l);
365
366 vector<size_t>* _patternLinks_father_son = &likelihoodData_->getArrayPositions(father->getId(), son->getId());
367 VVVdouble* _likelihoods_son = &likelihoodData_->getLikelihoodArray(son->getId());
368
369 if (son == branch)
370 {
371 VVVdouble* dpxy__son = &dpxy_[son->getId()];
372
373 for (size_t i = 0; i < nbSites; i++)
374 {
375 VVdouble* _likelihoods_son_i = &(*_likelihoods_son)[(*_patternLinks_father_son)[i]];
376 VVdouble* _dLikelihoods_father_i = &(*_dLikelihoods_father)[i];
377 for (size_t c = 0; c < nbClasses_; c++)
378 {
379 Vdouble* _likelihoods_son_i_c = &(*_likelihoods_son_i)[c];
380 Vdouble* _dLikelihoods_father_i_c = &(*_dLikelihoods_father_i)[c];
381 VVdouble* dpxy__son_c = &(*dpxy__son)[c];
382 for (size_t x = 0; x < nbStates_; x++)
383 {
384 double dl = 0;
385 Vdouble* dpxy__son_c_x = &(*dpxy__son_c)[x];
386 for (size_t y = 0; y < nbStates_; y++)
387 {
388 dl += (*dpxy__son_c_x)[y] * (*_likelihoods_son_i_c)[y];
389 }
390 (*_dLikelihoods_father_i_c)[x] *= dl;
391 }
392 }
393 }
394 }
395 else
396 {
397 VVVdouble* pxy__son = &pxy_[son->getId()];
398 for (size_t i = 0; i < nbSites; i++)
399 {
400 VVdouble* _likelihoods_son_i = &(*_likelihoods_son)[(*_patternLinks_father_son)[i]];
401 VVdouble* _dLikelihoods_father_i = &(*_dLikelihoods_father)[i];
402 for (size_t c = 0; c < nbClasses_; c++)
403 {
404 Vdouble* _likelihoods_son_i_c = &(*_likelihoods_son_i)[c];
405 Vdouble* _dLikelihoods_father_i_c = &(*_dLikelihoods_father_i)[c];
406 VVdouble* pxy__son_c = &(*pxy__son)[c];
407 for (size_t x = 0; x < nbStates_; x++)
408 {
409 double dl = 0;
410 Vdouble* pxy__son_c_x = &(*pxy__son_c)[x];
411 for (size_t y = 0; y < nbStates_; y++)
412 {
413 dl += (*pxy__son_c_x)[y] * (*_likelihoods_son_i_c)[y];
414 }
415 (*_dLikelihoods_father_i_c)[x] *= dl;
416 }
417 }
418 }
419 }
420 }
421
422 // Now we go down the tree toward the root node:
424}
425
426/******************************************************************************/
427
429{
430 const Node* father = node->getFather();
431 // We assume that the _dLikelihoods array has been filled for the current node 'node'.
432 // We will evaluate the array for the father node.
433 if (father == NULL)
434 return; // We reached the root!
435
436 // Compute dLikelihoods array for the father node.
437 // Fist initialize to 1:
438 VVVdouble* _dLikelihoods_father = &likelihoodData_->getDLikelihoodArray(father->getId());
439 size_t nbSites = _dLikelihoods_father->size();
440 for (size_t i = 0; i < nbSites; i++)
441 {
442 VVdouble* _dLikelihoods_father_i = &(*_dLikelihoods_father)[i];
443 for (size_t c = 0; c < nbClasses_; c++)
444 {
445 Vdouble* _dLikelihoods_father_i_c = &(*_dLikelihoods_father_i)[c];
446 for (size_t s = 0; s < nbStates_; s++)
447 {
448 (*_dLikelihoods_father_i_c)[s] = 1.;
449 }
450 }
451 }
452
453 size_t nbNodes = father->getNumberOfSons();
454 for (size_t l = 0; l < nbNodes; l++)
455 {
456 const Node* son = father->getSon(l);
457 VVVdouble* pxy__son = &pxy_[son->getId()];
458
459 vector<size_t>* _patternLinks_father_son = &likelihoodData_->getArrayPositions(father->getId(), son->getId());
460
461 if (son == node)
462 {
463 VVVdouble* _dLikelihoods_son = &likelihoodData_->getDLikelihoodArray(son->getId());
464 for (size_t i = 0; i < nbSites; i++)
465 {
466 VVdouble* _dLikelihoods_son_i = &(*_dLikelihoods_son)[(*_patternLinks_father_son)[i]];
467 VVdouble* _dLikelihoods_father_i = &(*_dLikelihoods_father)[i];
468 for (size_t c = 0; c < nbClasses_; c++)
469 {
470 Vdouble* _dLikelihoods_son_i_c = &(*_dLikelihoods_son_i)[c];
471 Vdouble* _dLikelihoods_father_i_c = &(*_dLikelihoods_father_i)[c];
472 VVdouble* pxy__son_c = &(*pxy__son)[c];
473 for (size_t x = 0; x < nbStates_; x++)
474 {
475 double dl = 0;
476 Vdouble* pxy__son_c_x = &(*pxy__son_c)[x];
477 for (size_t y = 0; y < nbStates_; y++)
478 {
479 dl += (*pxy__son_c_x)[y] * (*_dLikelihoods_son_i_c)[y];
480 }
481 (*_dLikelihoods_father_i_c)[x] *= dl;
482 }
483 }
484 }
485 }
486 else
487 {
488 VVVdouble* _likelihoods_son = &likelihoodData_->getLikelihoodArray(son->getId());
489 for (size_t i = 0; i < nbSites; i++)
490 {
491 VVdouble* _likelihoods_son_i = &(*_likelihoods_son)[(*_patternLinks_father_son)[i]];
492 VVdouble* _dLikelihoods_father_i = &(*_dLikelihoods_father)[i];
493 for (size_t c = 0; c < nbClasses_; c++)
494 {
495 Vdouble* _likelihoods_son_i_c = &(*_likelihoods_son_i)[c];
496 Vdouble* _dLikelihoods_father_i_c = &(*_dLikelihoods_father_i)[c];
497 VVdouble* pxy__son_c = &(*pxy__son)[c];
498 for (size_t x = 0; x < nbStates_; x++)
499 {
500 double dl = 0;
501 Vdouble* pxy__son_c_x = &(*pxy__son_c)[x];
502 for (size_t y = 0; y < nbStates_; y++)
503 {
504 dl += (*pxy__son_c_x)[y] * (*_likelihoods_son_i_c)[y];
505 }
506 (*_dLikelihoods_father_i_c)[x] *= dl;
507 }
508 }
509 }
510 }
511 }
512
513 // Next step: move toward grand father...
515}
516
517/******************************************************************************
518 * Second Order Derivatives *
519 ******************************************************************************/
521 size_t site,
522 size_t rateClass) const
523{
524 double d2l = 0;
525 Vdouble* d2la = &likelihoodData_->getD2LikelihoodArray(tree_->getRootNode()->getId())[likelihoodData_->getRootArrayPosition(site)][rateClass];
526 for (size_t i = 0; i < nbStates_; i++)
527 {
528 d2l += (*d2la)[i] * rootFreqs_[i];
529 }
530 return d2l;
531}
532
533/******************************************************************************/
534
536{
537 // Derivative of the sum is the sum of derivatives:
538 double d2l = 0;
539 for (size_t i = 0; i < nbClasses_; i++)
540 {
541 d2l += getD2LikelihoodForASiteForARateClass(site, i) * rateDistribution_->getProbability(i);
542 }
543 return d2l;
544}
545
546/******************************************************************************/
547
549{
552}
553
554/******************************************************************************/
555
557{
558 // Derivative of the sum is the sum of derivatives:
559 double dl = 0;
560 for (size_t i = 0; i < nbSites_; i++)
561 {
563 }
564 return dl;
565}
566
567/******************************************************************************/
568
570{
571 if (!hasParameter(variable))
572 throw ParameterNotFoundException("RHomogeneousTreeLikelihood::getSecondOrderDerivative().", variable);
574 {
575 throw Exception("Derivatives respective to rate distribution parameter are not implemented.");
576 }
578 {
579 throw Exception("Derivatives respective to substitution model parameters are not implemented.");
580 }
581
582 const_cast<RHomogeneousTreeLikelihood*>(this)->computeTreeD2Likelihood(variable);
583
584 return -getD2LogLikelihood();
585}
586
587/******************************************************************************/
588
590{
591 // Get the node with the branch whose length must be derivated:
592 size_t brI = TextTools::to<size_t>(variable.substr(5));
593 const Node* branch = nodes_[brI];
594 const Node* father = branch->getFather();
595
596 // Compute dLikelihoods array for the father node.
597 // Fist initialize to 1:
598 VVVdouble* _d2Likelihoods_father = &likelihoodData_->getD2LikelihoodArray(father->getId());
599 size_t nbSites = _d2Likelihoods_father->size();
600 for (size_t i = 0; i < nbSites; i++)
601 {
602 VVdouble* _d2Likelihoods_father_i = &(*_d2Likelihoods_father)[i];
603 for (size_t c = 0; c < nbClasses_; c++)
604 {
605 Vdouble* _d2Likelihoods_father_i_c = &(*_d2Likelihoods_father_i)[c];
606 for (size_t s = 0; s < nbStates_; s++)
607 {
608 (*_d2Likelihoods_father_i_c)[s] = 1.;
609 }
610 }
611 }
612
613 size_t nbNodes = father->getNumberOfSons();
614 for (size_t l = 0; l < nbNodes; l++)
615 {
616 const Node* son = father->getSon(l);
617
618 vector<size_t>* _patternLinks_father_son = &likelihoodData_->getArrayPositions(father->getId(), son->getId());
619 VVVdouble* _likelihoods_son = &likelihoodData_->getLikelihoodArray(son->getId());
620
621 if (son == branch)
622 {
623 VVVdouble* d2pxy__son = &d2pxy_[son->getId()];
624 for (size_t i = 0; i < nbSites; i++)
625 {
626 VVdouble* _likelihoods_son_i = &(*_likelihoods_son)[(*_patternLinks_father_son)[i]];
627 VVdouble* _d2Likelihoods_father_i = &(*_d2Likelihoods_father)[i];
628 for (size_t c = 0; c < nbClasses_; c++)
629 {
630 Vdouble* _likelihoods_son_i_c = &(*_likelihoods_son_i)[c];
631 Vdouble* _d2Likelihoods_father_i_c = &(*_d2Likelihoods_father_i)[c];
632 VVdouble* d2pxy__son_c = &(*d2pxy__son)[c];
633 for (size_t x = 0; x < nbStates_; x++)
634 {
635 double d2l = 0;
636 Vdouble* d2pxy__son_c_x = &(*d2pxy__son_c)[x];
637 for (size_t y = 0; y < nbStates_; y++)
638 {
639 d2l += (*d2pxy__son_c_x)[y] * (*_likelihoods_son_i_c)[y];
640 }
641 (*_d2Likelihoods_father_i_c)[x] *= d2l;
642 }
643 }
644 }
645 }
646 else
647 {
648 VVVdouble* pxy__son = &pxy_[son->getId()];
649 for (size_t i = 0; i < nbSites; i++)
650 {
651 VVdouble* _likelihoods_son_i = &(*_likelihoods_son)[(*_patternLinks_father_son)[i]];
652 VVdouble* _d2Likelihoods_father_i = &(*_d2Likelihoods_father)[i];
653 for (size_t c = 0; c < nbClasses_; c++)
654 {
655 Vdouble* _likelihoods_son_i_c = &(*_likelihoods_son_i)[c];
656 Vdouble* _d2Likelihoods_father_i_c = &(*_d2Likelihoods_father_i)[c];
657 VVdouble* pxy__son_c = &(*pxy__son)[c];
658 for (size_t x = 0; x < nbStates_; x++)
659 {
660 double d2l = 0;
661 Vdouble* pxy__son_c_x = &(*pxy__son_c)[x];
662 for (size_t y = 0; y < nbStates_; y++)
663 {
664 d2l += (*pxy__son_c_x)[y] * (*_likelihoods_son_i_c)[y];
665 }
666 (*_d2Likelihoods_father_i_c)[x] *= d2l;
667 }
668 }
669 }
670 }
671 }
672
673 // Now we go down the tree toward the root node:
675}
676
677/******************************************************************************/
678
680{
681 const Node* father = node->getFather();
682 // We assume that the _dLikelihoods array has been filled for the current node 'node'.
683 // We will evaluate the array for the father node.
684 if (father == NULL)
685 return; // We reached the root!
686
687 // Compute dLikelihoods array for the father node.
688 // Fist initialize to 1:
689 VVVdouble* _d2Likelihoods_father = &likelihoodData_->getD2LikelihoodArray(father->getId());
690 size_t nbSites = _d2Likelihoods_father->size();
691 for (size_t i = 0; i < nbSites; i++)
692 {
693 VVdouble* _d2Likelihoods_father_i = &(*_d2Likelihoods_father)[i];
694 for (size_t c = 0; c < nbClasses_; c++)
695 {
696 Vdouble* _d2Likelihoods_father_i_c = &(*_d2Likelihoods_father_i)[c];
697 for (size_t s = 0; s < nbStates_; s++)
698 {
699 (*_d2Likelihoods_father_i_c)[s] = 1.;
700 }
701 }
702 }
703
704 size_t nbNodes = father->getNumberOfSons();
705 for (size_t l = 0; l < nbNodes; l++)
706 {
707 const Node* son = father->getSon(l);
708
709 VVVdouble* pxy__son = &pxy_[son->getId()];
710 vector<size_t>* _patternLinks_father_son = &likelihoodData_->getArrayPositions(father->getId(), son->getId());
711
712 if (son == node)
713 {
714 VVVdouble* _d2Likelihoods_son = &likelihoodData_->getD2LikelihoodArray(son->getId());
715 for (size_t i = 0; i < nbSites; i++)
716 {
717 VVdouble* _d2Likelihoods_son_i = &(*_d2Likelihoods_son)[(*_patternLinks_father_son)[i]];
718 VVdouble* _d2Likelihoods_father_i = &(*_d2Likelihoods_father)[i];
719 for (size_t c = 0; c < nbClasses_; c++)
720 {
721 Vdouble* _d2Likelihoods_son_i_c = &(*_d2Likelihoods_son_i)[c];
722 Vdouble* _d2Likelihoods_father_i_c = &(*_d2Likelihoods_father_i)[c];
723 VVdouble* pxy__son_c = &(*pxy__son)[c];
724 for (size_t x = 0; x < nbStates_; x++)
725 {
726 double d2l = 0;
727 Vdouble* pxy__son_c_x = &(*pxy__son_c)[x];
728 for (size_t y = 0; y < nbStates_; y++)
729 {
730 d2l += (*pxy__son_c_x)[y] * (*_d2Likelihoods_son_i_c)[y];
731 }
732 (*_d2Likelihoods_father_i_c)[x] *= d2l;
733 }
734 }
735 }
736 }
737 else
738 {
739 VVVdouble* _likelihoods_son = &likelihoodData_->getLikelihoodArray(son->getId());
740 for (size_t i = 0; i < nbSites; i++)
741 {
742 VVdouble* _likelihoods_son_i = &(*_likelihoods_son)[(*_patternLinks_father_son)[i]];
743 VVdouble* _d2Likelihoods_father_i = &(*_d2Likelihoods_father)[i];
744 for (size_t c = 0; c < nbClasses_; c++)
745 {
746 Vdouble* _likelihoods_son_i_c = &(*_likelihoods_son_i)[c];
747 Vdouble* _d2Likelihoods_father_i_c = &(*_d2Likelihoods_father_i)[c];
748 VVdouble* pxy__son_c = &(*pxy__son)[c];
749 for (size_t x = 0; x < nbStates_; x++)
750 {
751 double dl = 0;
752 Vdouble* pxy__son_c_x = &(*pxy__son_c)[x];
753 for (size_t y = 0; y < nbStates_; y++)
754 {
755 dl += (*pxy__son_c_x)[y] * (*_likelihoods_son_i_c)[y];
756 }
757 (*_d2Likelihoods_father_i_c)[x] *= dl;
758 }
759 }
760 }
761 }
762 }
763
764 // Next step: move toward grand father...
766}
767
768/******************************************************************************/
769
771{
772 computeSubtreeLikelihood(tree_->getRootNode());
773}
774
775/******************************************************************************/
776
778{
779 if (node->isLeaf())
780 return;
781
782 size_t nbSites = likelihoodData_->getLikelihoodArray(node->getId()).size();
783 size_t nbNodes = node->getNumberOfSons();
784
785 // Must reset the likelihood array first (i.e. set all of them to 1):
786 VVVdouble* _likelihoods_node = &likelihoodData_->getLikelihoodArray(node->getId());
787 for (size_t i = 0; i < nbSites; i++)
788 {
789 // For each site in the sequence,
790 VVdouble* _likelihoods_node_i = &(*_likelihoods_node)[i];
791 for (size_t c = 0; c < nbClasses_; c++)
792 {
793 // For each rate class,
794 Vdouble* _likelihoods_node_i_c = &(*_likelihoods_node_i)[c];
795 for (size_t x = 0; x < nbStates_; x++)
796 {
797 // For each initial state,
798 (*_likelihoods_node_i_c)[x] = 1.;
799 }
800 }
801 }
802
803 for (size_t l = 0; l < nbNodes; l++)
804 {
805 // For each son node,
806
807 const Node* son = node->getSon(l);
808
809 computeSubtreeLikelihood(son); // Recursive method:
810
811 VVVdouble* pxy__son = &pxy_[son->getId()];
812 vector<size_t>* _patternLinks_node_son = &likelihoodData_->getArrayPositions(node->getId(), son->getId());
813 VVVdouble* _likelihoods_son = &likelihoodData_->getLikelihoodArray(son->getId());
814
815 for (size_t i = 0; i < nbSites; i++)
816 {
817 // For each site in the sequence,
818 VVdouble* _likelihoods_son_i = &(*_likelihoods_son)[(*_patternLinks_node_son)[i]];
819 VVdouble* _likelihoods_node_i = &(*_likelihoods_node)[i];
820 for (size_t c = 0; c < nbClasses_; c++)
821 {
822 // For each rate class,
823 Vdouble* _likelihoods_son_i_c = &(*_likelihoods_son_i)[c];
824 Vdouble* _likelihoods_node_i_c = &(*_likelihoods_node_i)[c];
825 VVdouble* pxy__son_c = &(*pxy__son)[c];
826 for (size_t x = 0; x < nbStates_; x++)
827 {
828 // For each initial state,
829 Vdouble* pxy__son_c_x = &(*pxy__son_c)[x];
830 double likelihood = 0;
831 for (size_t y = 0; y < nbStates_; y++)
832 {
833 likelihood += (*pxy__son_c_x)[y] * (*_likelihoods_son_i_c)[y];
834 }
835
836 (*_likelihoods_node_i_c)[x] *= likelihood;
837 }
838 }
839 }
840 }
841}
842
843/******************************************************************************/
844
846{
847 cout << "Likelihoods at node " << node->getName() << ": " << endl;
849 cout << " ***" << endl;
850}
851
852/*******************************************************************************/
static void displayLikelihoodArray(const VVVdouble &likelihoodArray)
Print the likelihood array to terminal (debugging tool).
Partial implementation for homogeneous model of the TreeLikelihood interface.
std::vector< Node * > nodes_
Pointer toward all nodes in the tree.
virtual void applyParameters()
All parameters are stored in a parameter list. This function apply these parameters to the substituti...
std::shared_ptr< TransitionModelInterface > model_
virtual void computeTransitionProbabilitiesForNode(const Node *node)
Fill the pxy_, dpxy_ and d2pxy_ arrays for one node.
ParameterList getRateDistributionParameters() const
Get the parameters associated to the rate distribution.
ParameterList getSubstitutionModelParameters() const
Get the parameters associated to substitution model(s).
virtual void computeAllTransitionProbabilities()
Fill the pxy_, dpxy_ and d2pxy_ arrays for all nodes.
AbstractHomogeneousTreeLikelihood & operator=(const AbstractHomogeneousTreeLikelihood &lik)
Assignation operator.
void setParametersValues(const ParameterList &parameters) override
bool hasParameter(const std::string &name) const override
const AlignmentDataInterface & data() const
Get the dataset for which the likelihood must be evaluated.
std::unique_ptr< const AlignmentDataInterface > data_
std::shared_ptr< TreeTemplate< Node > > tree_
static void displayTask(const std::string &text, bool eof=false)
static void displayTaskDone()
static void displayResult(const std::string &text, const T &result)
discrete Rate Across Sites, (simple) Recursive likelihood data structure.
void setTree(std::shared_ptr< const TreeTemplate< Node > > tree)
Set the tree associated to the data.
DRASRTreeLikelihoodData * clone() const
VVVdouble & getD2LikelihoodArray(int nodeId)
void initLikelihoods(const AlignmentDataInterface &sites, const TransitionModelInterface &model)
size_t getRootArrayPosition(size_t currentPosition) const
VVVdouble & getLikelihoodArray(int nodeId)
VVVdouble & getDLikelihoodArray(int nodeId)
const std::vector< size_t > & getArrayPositions(int parentId, int sonId) const
The phylogenetic node class.
Definition: Node.h:59
virtual std::string getName() const
Get the name associated to this node, if there is one, otherwise throw a NodeException.
Definition: Node.h:203
virtual int getId() const
Get the node's id.
Definition: Node.h:170
virtual const Node * getFather() const
Get the father of this node is there is one.
Definition: Node.h:306
virtual const Node * getSon(size_t pos) const
Definition: Node.h:362
virtual bool isLeaf() const
Definition: Node.h:679
virtual size_t getNumberOfSons() const
Definition: Node.h:355
size_t size() const
static std::unique_ptr< AlignmentDataInterface > getSequenceSubset(const AlignmentDataInterface &sequenceSet, const std::shared_ptr< N > node, const AssociationTreeGraphImplObserver< N, E, I > &tree)
Extract the sequences corresponding to a given subtree.
Definition: PatternTools.h:44
This class implement the 'traditional' way of computing likelihood for a tree.
virtual double getD2LikelihoodForASiteForARateClass(size_t site, size_t rateClass) const
virtual double getDLikelihoodForASiteForARateClass(size_t site, size_t rateClass) const
void setParameters(const ParameterList &parameters)
Implements the Function interface.
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.
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.
double getLogLikelihoodForASite(size_t site) const
Get the logarithm of the likelihood for a site.
RHomogeneousTreeLikelihood & operator=(const RHomogeneousTreeLikelihood &lik)
virtual void displayLikelihood(const Node *node)
This method is mainly for debugging purpose.
virtual void computeSubtreeLikelihood(const Node *node)
Compute the likelihood for a subtree defined by the Tree::Node node.
double getSecondOrderDerivative(const std::string &variable) const
double getLikelihoodForASite(size_t site) const
Get 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 setData(const AlignmentDataInterface &sites)
Set the dataset for which the likelihood must be evaluated.
virtual double getD2LikelihoodForASite(size_t site) const
virtual void computeTreeDLikelihood(const std::string &variable)
virtual void computeTreeD2Likelihood(const std::string &variable)
RHomogeneousTreeLikelihood(const Tree &tree, std::shared_ptr< TransitionModelInterface > model, std::shared_ptr< DiscreteDistributionInterface > rDist, bool checkRooted=true, bool verbose=true, bool usePatterns=true)
Build a new RHomogeneousTreeLikelihood object without data.
virtual double getDLikelihoodForASite(size_t site) const
double getLikelihood() const
Get the likelihood for the whole dataset.
double getFirstOrderDerivative(const std::string &variable) const
virtual double getD2LogLikelihoodForASite(size_t site) const
double getLikelihoodForASiteForARateClass(size_t site, size_t rateClass) const
Get the likelihood for a site knowing its rate class.
virtual double getDLogLikelihoodForASite(size_t site) const
double getLogLikelihood() const
Get the logarithm of the likelihood for the whole dataset.
void init_(bool usePatterns)
Method called by constructors.
DRASRTreeLikelihoodData * likelihoodData_
virtual void computeDownSubtreeD2Likelihood(const Node *)
void fireParameterChanged(const ParameterList &params)
virtual void computeDownSubtreeDLikelihood(const Node *)
Interface for phylogenetic tree objects.
Definition: Tree.h:115
std::string toString(T t)
Defines the basic types of data flow nodes.
double log(const ExtendedFloat &ef)
std::vector< double > Vdouble
ExtendedFloat pow(const ExtendedFloat &ef, double exp)
std::vector< VVdouble > VVVdouble
std::vector< Vdouble > VVdouble