bpp-phyl  2.2.0
RNonHomogeneousTreeLikelihood.cpp
Go to the documentation of this file.
1 //
2 // File: RNonHomogeneousTreeLikelihood.cpp
3 // Created by: Julien Dutheil
4 // Created on: Tue Oct 09 16:07 2007
5 // From file: RHomogeneousTreeLikelihood.cpp
6 //
7 
8 /*
9  Copyright or © or Copr. Bio++ Development Team, (November 16, 2004)
10 
11  This software is a computer program whose purpose is to provide classes
12  for phylogenetic data analysis.
13 
14  This software is governed by the CeCILL license under French law and
15  abiding by the rules of distribution of free software. You can use,
16  modify and/ or redistribute the software under the terms of the CeCILL
17  license as circulated by CEA, CNRS and INRIA at the following URL
18  "http://www.cecill.info".
19 
20  As a counterpart to the access to the source code and rights to copy,
21  modify and redistribute granted by the license, users are provided only
22  with a limited warranty and the software's author, the holder of the
23  economic rights, and the successive licensors have only limited
24  liability.
25 
26  In this respect, the user's attention is drawn to the risks associated
27  with loading, using, modifying and/or developing or reproducing the
28  software by the user in light of its specific status of free software,
29  that may mean that it is complicated to manipulate, and that also
30  therefore means that it is reserved for developers and experienced
31  professionals having in-depth computer knowledge. Users are therefore
32  encouraged to load and test the software's suitability as regards their
33  requirements in conditions enabling the security of their systems and/or
34  data to be ensured and, more generally, to use and operate it in the
35  same conditions as regards security.
36 
37  The fact that you are presently reading this means that you have had
38  knowledge of the CeCILL license and that you accept its terms.
39  */
40 
42 #include "../PatternTools.h"
43 
44 #include <Bpp/Text/TextTools.h>
45 #include <Bpp/App/ApplicationTools.h>
46 
47 using namespace bpp;
48 
49 // From the STL:
50 #include <iostream>
51 
52 using namespace std;
53 
54 /******************************************************************************/
55 
57  const Tree& tree,
58  SubstitutionModelSet* modelSet,
59  DiscreteDistribution* rDist,
60  bool verbose,
61  bool usePatterns,
62  bool reparametrizeRoot)
63 throw (Exception) :
64  AbstractNonHomogeneousTreeLikelihood(tree, modelSet, rDist, verbose, reparametrizeRoot),
65  likelihoodData_(0),
66  minusLogLik_(-1.)
67 {
68  if (!modelSet->isFullySetUpFor(tree))
69  throw Exception("RNonHomogeneousTreeLikelihood(constructor). Model set is not fully specified.");
70  init_(usePatterns);
71 }
72 
73 /******************************************************************************/
74 
76  const Tree& tree,
77  const SiteContainer& data,
78  SubstitutionModelSet* modelSet,
79  DiscreteDistribution* rDist,
80  bool verbose,
81  bool usePatterns,
82  bool reparametrizeRoot)
83 throw (Exception) :
84  AbstractNonHomogeneousTreeLikelihood(tree, modelSet, rDist, verbose, reparametrizeRoot),
85  likelihoodData_(0),
86  minusLogLik_(-1.)
87 {
88  if (!modelSet->isFullySetUpFor(tree))
89  throw Exception("RNonHomogeneousTreeLikelihood(constructor). Model set is not fully specified.");
90  init_(usePatterns);
91  setData(data);
92 }
93 
94 /******************************************************************************/
95 
96 void RNonHomogeneousTreeLikelihood::init_(bool usePatterns) throw (Exception)
97 {
98  likelihoodData_ = new DRASRTreeLikelihoodData(
99  tree_,
100  rateDistribution_->getNumberOfCategories(),
101  usePatterns);
102 }
103 
104 /******************************************************************************/
105 
107  const RNonHomogeneousTreeLikelihood& lik) :
109  likelihoodData_(0),
110  minusLogLik_(lik.minusLogLik_)
111 {
114 }
115 
116 /******************************************************************************/
117 
120 {
122  if (likelihoodData_) delete likelihoodData_;
126  return *this;
127 }
128 
129 /******************************************************************************/
130 
132 {
133  delete likelihoodData_;
134 }
135 
136 /******************************************************************************/
137 
138 void RNonHomogeneousTreeLikelihood::setData(const SiteContainer& sites) throw (Exception)
139 {
140  if (data_) delete data_;
141  data_ = PatternTools::getSequenceSubset(sites, *tree_->getRootNode());
142  if (verbose_) ApplicationTools::displayTask("Initializing data structure");
143  likelihoodData_->initLikelihoods(*data_, *modelSet_->getModel(0)); //We assume here that all models have the same number of states, and that they have the same 'init' method,
144  //Which is a reasonable assumption as long as they share the same alphabet.
145  if (verbose_) ApplicationTools::displayTaskDone();
146 
147  nbSites_ = likelihoodData_->getNumberOfSites();
148  nbDistinctSites_ = likelihoodData_->getNumberOfDistinctSites();
149  nbStates_ = likelihoodData_->getNumberOfStates();
150 
151  if (verbose_) ApplicationTools::displayResult("Number of distinct sites",
152  TextTools::toString(nbDistinctSites_));
153  initialized_ = false;
154 }
155 
156 /******************************************************************************/
157 
159 {
160  double l = 1.;
161  for (size_t i = 0; i < nbSites_; i++)
162  {
163  l *= getLikelihoodForASite(i);
164  }
165  return l;
166 }
167 
168 /******************************************************************************/
169 
171 {
172  double ll = 0;
173  vector<double> la(nbSites_);
174  for (size_t i = 0; i < nbSites_; i++)
175  {
176  la[i] = getLogLikelihoodForASite(i);
177  }
178  sort(la.begin(), la.end());
179  for (size_t i = nbSites_; i > 0; i--)
180  {
181  ll += la[i - 1];
182  }
183  return ll;
184 }
185 
186 /******************************************************************************/
187 
189 {
190  double l = 0;
191  for (size_t i = 0; i < nbClasses_; i++)
192  {
193  l += getLikelihoodForASiteForARateClass(site, i) * rateDistribution_->getProbability(i);
194  }
195  return l;
196 }
197 
198 /******************************************************************************/
199 
201 {
202  double l = 0;
203  for (size_t i = 0; i < nbClasses_; i++)
204  {
205  l += getLikelihoodForASiteForARateClass(site, i) * rateDistribution_->getProbability(i);
206  }
207  //if(l <= 0.) cerr << "WARNING!!! Negative likelihood." << endl;
208  if (l < 0) l = 0; //May happen because of numerical errors.
209  return log(l);
210 }
211 
212 /******************************************************************************/
213 
214 double RNonHomogeneousTreeLikelihood::getLikelihoodForASiteForARateClass(size_t site, size_t rateClass) const
215 {
216  double l = 0;
217  Vdouble* la = &likelihoodData_->getLikelihoodArray(tree_->getRootNode()->getId())[likelihoodData_->getRootArrayPosition(site)][rateClass];
218  for (size_t i = 0; i < nbStates_; i++)
219  {
220  l += (*la)[i] * rootFreqs_[i];
221  }
222  return l;
223 }
224 
225 /******************************************************************************/
226 
228 {
229  double l = 0;
230  Vdouble* la = &likelihoodData_->getLikelihoodArray(tree_->getRootNode()->getId())[likelihoodData_->getRootArrayPosition(site)][rateClass];
231  for (size_t i = 0; i < nbStates_; i++)
232  {
233  l += (*la)[i] * rootFreqs_[i];
234  }
235  return log(l);
236 }
237 
238 /******************************************************************************/
239 
240 double RNonHomogeneousTreeLikelihood::getLikelihoodForASiteForARateClassForAState(size_t site, size_t rateClass, int state) const
241 {
242  return likelihoodData_->getLikelihoodArray(tree_->getRootNode()->getId())[likelihoodData_->getRootArrayPosition(site)][rateClass][static_cast<size_t>(state)];
243 }
244 
245 /******************************************************************************/
246 
247 double RNonHomogeneousTreeLikelihood::getLogLikelihoodForASiteForARateClassForAState(size_t site, size_t rateClass, int state) const
248 {
249  return log(likelihoodData_->getLikelihoodArray(tree_->getRootNode()->getId())[likelihoodData_->getRootArrayPosition(site)][rateClass][static_cast<size_t>(state)]);
250 }
251 
252 /******************************************************************************/
253 
254 void RNonHomogeneousTreeLikelihood::setParameters(const ParameterList& parameters)
255 throw (ParameterNotFoundException, ConstraintException)
256 {
257  setParametersValues(parameters);
258 }
259 
260 /******************************************************************************/
261 
263 {
264  applyParameters();
265 
266  if (params.getCommonParametersWith(rateDistribution_->getIndependentParameters()).size() > 0)
267  {
269  }
270  else
271  {
272  vector<int> ids;
273  vector<string> tmp = params.getCommonParametersWith(modelSet_->getNodeParameters()).getParameterNames();
274  for (size_t i = 0; i < tmp.size(); i++)
275  {
276  vector<int> tmpv = modelSet_->getNodesWithParameter(tmp[i]);
277  ids = VectorTools::vectorUnion(ids, tmpv);
278  }
279  tmp = params.getCommonParametersWith(brLenParameters_).getParameterNames();
280  vector<const Node*> nodes;
281  for (size_t i = 0; i < ids.size(); i++)
282  {
283  nodes.push_back(idToNode_[ids[i]]);
284  }
285  vector<const Node*> tmpv;
286  bool test = false;
287  for (size_t i = 0; i < tmp.size(); i++)
288  {
289  if (tmp[i] == "BrLenRoot" || tmp[i] == "RootPosition")
290  {
291  if (!test)
292  {
293  tmpv.push_back(tree_->getRootNode()->getSon(0));
294  tmpv.push_back(tree_->getRootNode()->getSon(1));
295  test = true; //Add only once.
296  }
297  }
298  else
299  tmpv.push_back(nodes_[TextTools::to < size_t > (tmp[i].substr(5))]);
300  }
301  nodes = VectorTools::vectorUnion(nodes, tmpv);
302 
303  for (size_t i = 0; i < nodes.size(); i++)
304  {
306  }
308  }
310 
312 }
313 
314 /******************************************************************************/
315 
317 throw (Exception)
318 {
319  if (!isInitialized()) throw Exception("RNonHomogeneousTreeLikelihood::getValue(). Instance is not initialized.");
320  return minusLogLik_;
321 }
322 
323 /******************************************************************************
324 * First Order Derivatives *
325 ******************************************************************************/
326 
328  size_t site,
329  size_t rateClass) const
330 {
331  double dl = 0;
332  Vdouble* dla = &likelihoodData_->getDLikelihoodArray(tree_->getRootNode()->getId())[likelihoodData_->getRootArrayPosition(site)][rateClass];
333  for (size_t i = 0; i < nbStates_; i++)
334  {
335  dl += (*dla)[i] * rootFreqs_[i];
336  }
337  return dl;
338 }
339 
340 /******************************************************************************/
341 
343 {
344  // Derivative of the sum is the sum of derivatives:
345  double dl = 0;
346  for (size_t i = 0; i < nbClasses_; i++)
347  {
348  dl += getDLikelihoodForASiteForARateClass(site, i) * rateDistribution_->getProbability(i);
349  }
350  return dl;
351 }
352 
353 /******************************************************************************/
354 
356 {
357  // d(f(g(x)))/dx = dg(x)/dx . df(g(x))/dg :
358  return getDLikelihoodForASite(site) / getLikelihoodForASite(site);
359 }
360 
361 /******************************************************************************/
362 
364 {
365  // Derivative of the sum is the sum of derivatives:
366  double dl = 0;
367  for (size_t i = 0; i < nbSites_; i++)
368  {
369  dl += getDLogLikelihoodForASite(i);
370  }
371  return dl;
372 }
373 
374 /******************************************************************************/
375 
376 double RNonHomogeneousTreeLikelihood::getFirstOrderDerivative(const string& variable) const
377 throw (Exception)
378 {
379  if (!hasParameter(variable))
380  throw ParameterNotFoundException("RNonHomogeneousTreeLikelihood::getFirstOrderDerivative().", variable);
381  if (getRateDistributionParameters().hasParameter(variable))
382  {
383  throw Exception("Derivatives respective to rate distribution parameter are not implemented.");
384  }
385  if (getSubstitutionModelParameters().hasParameter(variable))
386  {
387  throw Exception("Derivatives respective to substitution model parameters are not implemented.");
388  }
389 
390  const_cast<RNonHomogeneousTreeLikelihood*>(this)->computeTreeDLikelihood(variable);
391  return -getDLogLikelihood();
392 }
393 
394 /******************************************************************************/
395 
397 {
398  if (variable == "BrLenRoot")
399  {
400  const Node* father = tree_->getRootNode();
401 
402  // Compute dLikelihoods array for the father node.
403  // Fist initialize to 1:
404  VVVdouble* _dLikelihoods_father = &likelihoodData_->getDLikelihoodArray(father->getId());
405  size_t nbSites = _dLikelihoods_father->size();
406  for (size_t i = 0; i < nbSites; i++)
407  {
408  VVdouble* _dLikelihoods_father_i = &(*_dLikelihoods_father)[i];
409  for (size_t c = 0; c < nbClasses_; c++)
410  {
411  Vdouble* _dLikelihoods_father_i_c = &(*_dLikelihoods_father_i)[c];
412  for (size_t s = 0; s < nbStates_; s++)
413  {
414  (*_dLikelihoods_father_i_c)[s] = 1.;
415  }
416  }
417  }
418 
419  size_t nbNodes = father->getNumberOfSons();
420  for (size_t l = 0; l < nbNodes; l++)
421  {
422  const Node* son = father->getSon(l);
423 
424  if (son->getId() == root1_)
425  {
426  const Node* root1 = father->getSon(0);
427  const Node* root2 = father->getSon(1);
428  vector<size_t> * _patternLinks_fatherroot1_ = &likelihoodData_->getArrayPositions(father->getId(), root1->getId());
429  vector<size_t> * _patternLinks_fatherroot2_ = &likelihoodData_->getArrayPositions(father->getId(), root2->getId());
430  VVVdouble* _likelihoodsroot1_ = &likelihoodData_->getLikelihoodArray(root1->getId());
431  VVVdouble* _likelihoodsroot2_ = &likelihoodData_->getLikelihoodArray(root2->getId());
432  double pos = getParameterValue("RootPosition");
433 
434  VVVdouble* dpxy_root1_ = &dpxy_[root1_];
435  VVVdouble* dpxy_root2_ = &dpxy_[root2_];
436  VVVdouble* pxy_root1_ = &pxy_[root1_];
437  VVVdouble* pxy_root2_ = &pxy_[root2_];
438  for (size_t i = 0; i < nbSites; i++)
439  {
440  VVdouble* _likelihoodsroot1__i = &(*_likelihoodsroot1_)[(*_patternLinks_fatherroot1_)[i]];
441  VVdouble* _likelihoodsroot2__i = &(*_likelihoodsroot2_)[(*_patternLinks_fatherroot2_)[i]];
442  VVdouble* _dLikelihoods_father_i = &(*_dLikelihoods_father)[i];
443  for (size_t c = 0; c < nbClasses_; c++)
444  {
445  Vdouble* _likelihoodsroot1__i_c = &(*_likelihoodsroot1__i)[c];
446  Vdouble* _likelihoodsroot2__i_c = &(*_likelihoodsroot2__i)[c];
447  Vdouble* _dLikelihoods_father_i_c = &(*_dLikelihoods_father_i)[c];
448  VVdouble* dpxy_root1__c = &(*dpxy_root1_)[c];
449  VVdouble* dpxy_root2__c = &(*dpxy_root2_)[c];
450  VVdouble* pxy_root1__c = &(*pxy_root1_)[c];
451  VVdouble* pxy_root2__c = &(*pxy_root2_)[c];
452  for (size_t x = 0; x < nbStates_; x++)
453  {
454  Vdouble* dpxy_root1__c_x = &(*dpxy_root1__c)[x];
455  Vdouble* dpxy_root2__c_x = &(*dpxy_root2__c)[x];
456  Vdouble* pxy_root1__c_x = &(*pxy_root1__c)[x];
457  Vdouble* pxy_root2__c_x = &(*pxy_root2__c)[x];
458  double dl1 = 0, dl2 = 0, l1 = 0, l2 = 0;
459  for (size_t y = 0; y < nbStates_; y++)
460  {
461  dl1 += (*dpxy_root1__c_x)[y] * (*_likelihoodsroot1__i_c)[y];
462  dl2 += (*dpxy_root2__c_x)[y] * (*_likelihoodsroot2__i_c)[y];
463  l1 += (*pxy_root1__c_x)[y] * (*_likelihoodsroot1__i_c)[y];
464  l2 += (*pxy_root2__c_x)[y] * (*_likelihoodsroot2__i_c)[y];
465  }
466  double dl = pos * dl1 * l2 + (1. - pos) * dl2 * l1;
467  (*_dLikelihoods_father_i_c)[x] *= dl;
468  }
469  }
470  }
471  }
472  else if (son->getId() == root2_)
473  {
474  //Do nothing, this was accounted when dealing with root1_
475  }
476  else
477  {
478  //Account for a putative multifurcation:
479  vector<size_t> * _patternLinks_father_son = &likelihoodData_->getArrayPositions(father->getId(), son->getId());
480  VVVdouble* _likelihoods_son = &likelihoodData_->getLikelihoodArray(son->getId());
481 
482  VVVdouble* pxy__son = &pxy_[son->getId()];
483  for (size_t i = 0; i < nbSites; i++)
484  {
485  VVdouble* _likelihoods_son_i = &(*_likelihoods_son)[(*_patternLinks_father_son)[i]];
486  VVdouble* _dLikelihoods_father_i = &(*_dLikelihoods_father)[i];
487  for (size_t c = 0; c < nbClasses_; c++)
488  {
489  Vdouble* _likelihoods_son_i_c = &(*_likelihoods_son_i)[c];
490  Vdouble* _dLikelihoods_father_i_c = &(*_dLikelihoods_father_i)[c];
491  VVdouble* pxy__son_c = &(*pxy__son)[c];
492  for (size_t x = 0; x < nbStates_; x++)
493  {
494  double dl = 0;
495  Vdouble* pxy__son_c_x = &(*pxy__son_c)[x];
496  for (size_t y = 0; y < nbStates_; y++)
497  {
498  dl += (*pxy__son_c_x)[y] * (*_likelihoods_son_i_c)[y];
499  }
500  (*_dLikelihoods_father_i_c)[x] *= dl;
501  }
502  }
503  }
504  }
505  }
506  return;
507  }
508  else if (variable == "RootPosition")
509  {
510  const Node* father = tree_->getRootNode();
511 
512  // Compute dLikelihoods array for the father node.
513  // Fist initialize to 1:
514  VVVdouble* _dLikelihoods_father = &likelihoodData_->getDLikelihoodArray(father->getId());
515  size_t nbSites = _dLikelihoods_father->size();
516  for (size_t i = 0; i < nbSites; i++)
517  {
518  VVdouble* _dLikelihoods_father_i = &(*_dLikelihoods_father)[i];
519  for (size_t c = 0; c < nbClasses_; c++)
520  {
521  Vdouble* _dLikelihoods_father_i_c = &(*_dLikelihoods_father_i)[c];
522  for (size_t s = 0; s < nbStates_; s++)
523  {
524  (*_dLikelihoods_father_i_c)[s] = 1.;
525  }
526  }
527  }
528 
529  size_t nbNodes = father->getNumberOfSons();
530  for (size_t l = 0; l < nbNodes; l++)
531  {
532  const Node* son = father->getSon(l);
533 
534  if (son->getId() == root1_)
535  {
536  const Node* root1 = father->getSon(0);
537  const Node* root2 = father->getSon(1);
538  vector<size_t> * _patternLinks_fatherroot1_ = &likelihoodData_->getArrayPositions(father->getId(), root1->getId());
539  vector<size_t> * _patternLinks_fatherroot2_ = &likelihoodData_->getArrayPositions(father->getId(), root2->getId());
540  VVVdouble* _likelihoodsroot1_ = &likelihoodData_->getLikelihoodArray(root1->getId());
541  VVVdouble* _likelihoodsroot2_ = &likelihoodData_->getLikelihoodArray(root2->getId());
542  double len = getParameterValue("BrLenRoot");
543 
544  VVVdouble* dpxy_root1_ = &dpxy_[root1_];
545  VVVdouble* dpxy_root2_ = &dpxy_[root2_];
546  VVVdouble* pxy_root1_ = &pxy_[root1_];
547  VVVdouble* pxy_root2_ = &pxy_[root2_];
548  for (size_t i = 0; i < nbSites; i++)
549  {
550  VVdouble* _likelihoodsroot1__i = &(*_likelihoodsroot1_)[(*_patternLinks_fatherroot1_)[i]];
551  VVdouble* _likelihoodsroot2__i = &(*_likelihoodsroot2_)[(*_patternLinks_fatherroot2_)[i]];
552  VVdouble* _dLikelihoods_father_i = &(*_dLikelihoods_father)[i];
553  for (size_t c = 0; c < nbClasses_; c++)
554  {
555  Vdouble* _likelihoodsroot1__i_c = &(*_likelihoodsroot1__i)[c];
556  Vdouble* _likelihoodsroot2__i_c = &(*_likelihoodsroot2__i)[c];
557  Vdouble* _dLikelihoods_father_i_c = &(*_dLikelihoods_father_i)[c];
558  VVdouble* dpxy_root1__c = &(*dpxy_root1_)[c];
559  VVdouble* dpxy_root2__c = &(*dpxy_root2_)[c];
560  VVdouble* pxy_root1__c = &(*pxy_root1_)[c];
561  VVdouble* pxy_root2__c = &(*pxy_root2_)[c];
562  for (size_t x = 0; x < nbStates_; x++)
563  {
564  Vdouble* dpxy_root1__c_x = &(*dpxy_root1__c)[x];
565  Vdouble* dpxy_root2__c_x = &(*dpxy_root2__c)[x];
566  Vdouble* pxy_root1__c_x = &(*pxy_root1__c)[x];
567  Vdouble* pxy_root2__c_x = &(*pxy_root2__c)[x];
568  double dl1 = 0, dl2 = 0, l1 = 0, l2 = 0;
569  for (size_t y = 0; y < nbStates_; y++)
570  {
571  dl1 += (*dpxy_root1__c_x)[y] * (*_likelihoodsroot1__i_c)[y];
572  dl2 += (*dpxy_root2__c_x)[y] * (*_likelihoodsroot2__i_c)[y];
573  l1 += (*pxy_root1__c_x)[y] * (*_likelihoodsroot1__i_c)[y];
574  l2 += (*pxy_root2__c_x)[y] * (*_likelihoodsroot2__i_c)[y];
575  }
576  double dl = len * (dl1 * l2 - dl2 * l1);
577  (*_dLikelihoods_father_i_c)[x] *= dl;
578  }
579  }
580  }
581  }
582  else if (son->getId() == root2_)
583  {
584  //Do nothing, this was accounted when dealing with root1_
585  }
586  else
587  {
588  //Account for a putative multifurcation:
589  vector<size_t> * _patternLinks_father_son = &likelihoodData_->getArrayPositions(father->getId(), son->getId());
590  VVVdouble* _likelihoods_son = &likelihoodData_->getLikelihoodArray(son->getId());
591 
592  VVVdouble* pxy__son = &pxy_[son->getId()];
593  for (size_t i = 0; i < nbSites; i++)
594  {
595  VVdouble* _likelihoods_son_i = &(*_likelihoods_son)[(*_patternLinks_father_son)[i]];
596  VVdouble* _dLikelihoods_father_i = &(*_dLikelihoods_father)[i];
597  for (size_t c = 0; c < nbClasses_; c++)
598  {
599  Vdouble* _likelihoods_son_i_c = &(*_likelihoods_son_i)[c];
600  Vdouble* _dLikelihoods_father_i_c = &(*_dLikelihoods_father_i)[c];
601  VVdouble* pxy__son_c = &(*pxy__son)[c];
602  for (size_t x = 0; x < nbStates_; x++)
603  {
604  double dl = 0;
605  Vdouble* pxy__son_c_x = &(*pxy__son_c)[x];
606  for (size_t y = 0; y < nbStates_; y++)
607  {
608  dl += (*pxy__son_c_x)[y] * (*_likelihoods_son_i_c)[y];
609  }
610  (*_dLikelihoods_father_i_c)[x] *= dl;
611  }
612  }
613  }
614  }
615  }
616  return;
617  }
618 
619  // Get the node with the branch whose length must be derivated:
620  size_t brI = TextTools::to<size_t>(variable.substr(5));
621  const Node* branch = nodes_[brI];
622  const Node* father = branch->getFather();
623  VVVdouble* _dLikelihoods_father = &likelihoodData_->getDLikelihoodArray(father->getId());
624 
625  // Compute dLikelihoods array for the father node.
626  // Fist initialize to 1:
627  size_t nbSites = _dLikelihoods_father->size();
628  for (size_t i = 0; i < nbSites; i++)
629  {
630  VVdouble* _dLikelihoods_father_i = &(*_dLikelihoods_father)[i];
631  for (size_t c = 0; c < nbClasses_; c++)
632  {
633  Vdouble* _dLikelihoods_father_i_c = &(*_dLikelihoods_father_i)[c];
634  for (size_t s = 0; s < nbStates_; s++)
635  {
636  (*_dLikelihoods_father_i_c)[s] = 1.;
637  }
638  }
639  }
640 
641  size_t nbNodes = father->getNumberOfSons();
642  for (size_t l = 0; l < nbNodes; l++)
643  {
644  const Node* son = father->getSon(l);
645 
646  vector<size_t> * _patternLinks_father_son = &likelihoodData_->getArrayPositions(father->getId(), son->getId());
647  VVVdouble* _likelihoods_son = &likelihoodData_->getLikelihoodArray(son->getId());
648 
649  if (son == branch)
650  {
651  VVVdouble* dpxy__son = &dpxy_[son->getId()];
652  for (size_t i = 0; i < nbSites; i++)
653  {
654  VVdouble* _likelihoods_son_i = &(*_likelihoods_son)[(*_patternLinks_father_son)[i]];
655  VVdouble* _dLikelihoods_father_i = &(*_dLikelihoods_father)[i];
656  for (size_t c = 0; c < nbClasses_; c++)
657  {
658  Vdouble* _likelihoods_son_i_c = &(*_likelihoods_son_i)[c];
659  Vdouble* _dLikelihoods_father_i_c = &(*_dLikelihoods_father_i)[c];
660  VVdouble* dpxy__son_c = &(*dpxy__son)[c];
661  for (size_t x = 0; x < nbStates_; x++)
662  {
663  double dl = 0;
664  Vdouble* dpxy__son_c_x = &(*dpxy__son_c)[x];
665  for (size_t y = 0; y < nbStates_; y++)
666  {
667  dl += (*dpxy__son_c_x)[y] * (*_likelihoods_son_i_c)[y];
668  }
669  (*_dLikelihoods_father_i_c)[x] *= dl;
670  }
671  }
672  }
673  }
674  else
675  {
676  VVVdouble* pxy__son = &pxy_[son->getId()];
677  for (size_t i = 0; i < nbSites; i++)
678  {
679  VVdouble* _likelihoods_son_i = &(*_likelihoods_son)[(*_patternLinks_father_son)[i]];
680  VVdouble* _dLikelihoods_father_i = &(*_dLikelihoods_father)[i];
681  for (size_t c = 0; c < nbClasses_; c++)
682  {
683  Vdouble* _likelihoods_son_i_c = &(*_likelihoods_son_i)[c];
684  Vdouble* _dLikelihoods_father_i_c = &(*_dLikelihoods_father_i)[c];
685  VVdouble* pxy__son_c = &(*pxy__son)[c];
686  for (size_t x = 0; x < nbStates_; x++)
687  {
688  double dl = 0;
689  Vdouble* pxy__son_c_x = &(*pxy__son_c)[x];
690  for (size_t y = 0; y < nbStates_; y++)
691  {
692  dl += (*pxy__son_c_x)[y] * (*_likelihoods_son_i_c)[y];
693  }
694  (*_dLikelihoods_father_i_c)[x] *= dl;
695  }
696  }
697  }
698  }
699  }
700 
701  // Now we go down the tree toward the root node:
703 }
704 
705 /******************************************************************************/
706 
708 {
709  const Node* father = node->getFather();
710  // We assume that the _dLikelihoods array has been filled for the current node 'node'.
711  // We will evaluate the array for the father node.
712  if (father == 0) return; // We reached the root!
713 
714  // Compute dLikelihoods array for the father node.
715  // Fist initialize to 1:
716  VVVdouble* _dLikelihoods_father = &likelihoodData_->getDLikelihoodArray(father->getId());
717  size_t nbSites = _dLikelihoods_father->size();
718  for (size_t i = 0; i < nbSites; i++)
719  {
720  VVdouble* _dLikelihoods_father_i = &(*_dLikelihoods_father)[i];
721  for (size_t c = 0; c < nbClasses_; c++)
722  {
723  Vdouble* _dLikelihoods_father_i_c = &(*_dLikelihoods_father_i)[c];
724  for (size_t s = 0; s < nbStates_; s++)
725  {
726  (*_dLikelihoods_father_i_c)[s] = 1.;
727  }
728  }
729  }
730 
731  size_t nbNodes = father->getNumberOfSons();
732  for (size_t l = 0; l < nbNodes; l++)
733  {
734  const Node* son = father->getSon(l);
735 
736  VVVdouble* pxy__son = &pxy_[son->getId()];
737  vector<size_t> * _patternLinks_father_son = &likelihoodData_->getArrayPositions(father->getId(), son->getId());
738 
739  if (son == node)
740  {
741  VVVdouble* _dLikelihoods_son = &likelihoodData_->getDLikelihoodArray(son->getId());
742  for (size_t i = 0; i < nbSites; i++)
743  {
744  VVdouble* _dLikelihoods_son_i = &(*_dLikelihoods_son)[(*_patternLinks_father_son)[i]];
745  VVdouble* _dLikelihoods_father_i = &(*_dLikelihoods_father)[i];
746  for (size_t c = 0; c < nbClasses_; c++)
747  {
748  Vdouble* _dLikelihoods_son_i_c = &(*_dLikelihoods_son_i)[c];
749  Vdouble* _dLikelihoods_father_i_c = &(*_dLikelihoods_father_i)[c];
750  VVdouble* pxy__son_c = &(*pxy__son)[c];
751  for (size_t x = 0; x < nbStates_; x++)
752  {
753  double dl = 0;
754  Vdouble* pxy__son_c_x = &(*pxy__son_c)[x];
755  for (size_t y = 0; y < nbStates_; y++)
756  {
757  dl += (*pxy__son_c_x)[y] * (*_dLikelihoods_son_i_c)[y];
758  }
759  (*_dLikelihoods_father_i_c)[x] *= dl;
760  }
761  }
762  }
763  }
764  else
765  {
766  VVVdouble* _likelihoods_son = &likelihoodData_->getLikelihoodArray(son->getId());
767  for (size_t i = 0; i < nbSites; i++)
768  {
769  VVdouble* _likelihoods_son_i = &(*_likelihoods_son)[(*_patternLinks_father_son)[i]];
770  VVdouble* _dLikelihoods_father_i = &(*_dLikelihoods_father)[i];
771  for (size_t c = 0; c < nbClasses_; c++)
772  {
773  Vdouble* _likelihoods_son_i_c = &(*_likelihoods_son_i)[c];
774  Vdouble* _dLikelihoods_father_i_c = &(*_dLikelihoods_father_i)[c];
775  VVdouble* pxy__son_c = &(*pxy__son)[c];
776  for (size_t x = 0; x < nbStates_; x++)
777  {
778  double dl = 0;
779  Vdouble* pxy__son_c_x = &(*pxy__son_c)[x];
780  for (size_t y = 0; y < nbStates_; y++)
781  {
782  dl += (*pxy__son_c_x)[y] * (*_likelihoods_son_i_c)[y];
783  }
784  (*_dLikelihoods_father_i_c)[x] *= dl;
785  }
786  }
787  }
788  }
789  }
790 
791  //Next step: move toward grand father...
793 }
794 
795 /******************************************************************************
796 * Second Order Derivatives *
797 ******************************************************************************/
798 
800  size_t site,
801  size_t rateClass) const
802 {
803  double d2l = 0;
804  Vdouble* d2la = &likelihoodData_->getD2LikelihoodArray(tree_->getRootNode()->getId())[likelihoodData_->getRootArrayPosition(site)][rateClass];
805  for (size_t i = 0; i < nbStates_; i++)
806  {
807  d2l += (*d2la)[i] * rootFreqs_[i];
808  }
809  return d2l;
810 }
811 
812 /******************************************************************************/
813 
815 {
816  // Derivative of the sum is the sum of derivatives:
817  double d2l = 0;
818  for (size_t i = 0; i < nbClasses_; i++)
819  {
820  d2l += getD2LikelihoodForASiteForARateClass(site, i) * rateDistribution_->getProbability(i);
821  }
822  return d2l;
823 }
824 
825 /******************************************************************************/
826 
828 {
830  - pow( getDLikelihoodForASite(site) / getLikelihoodForASite(site), 2);
831 }
832 
833 /******************************************************************************/
834 
836 {
837  // Derivative of the sum is the sum of derivatives:
838  double dl = 0;
839  for (size_t i = 0; i < nbSites_; i++)
840  {
842  }
843  return dl;
844 }
845 
846 /******************************************************************************/
847 
849 throw (Exception)
850 {
851  if (!hasParameter(variable))
852  throw ParameterNotFoundException("RNonHomogeneousTreeLikelihood::getSecondOrderDerivative().", variable);
853  if (getRateDistributionParameters().hasParameter(variable))
854  {
855  throw Exception("Derivatives respective to rate distribution parameter are not implemented.");
856  }
857  if (getSubstitutionModelParameters().hasParameter(variable))
858  {
859  throw Exception("Derivatives respective to substitution model parameters are not implemented.");
860  }
861 
862  const_cast<RNonHomogeneousTreeLikelihood*>(this)->computeTreeD2Likelihood(variable);
863  return -getD2LogLikelihood();
864 }
865 
866 /******************************************************************************/
867 
869 {
870  if (variable == "BrLenRoot")
871  {
872  const Node* father = tree_->getRootNode();
873 
874  // Compute dLikelihoods array for the father node.
875  // Fist initialize to 1:
876  VVVdouble* _d2Likelihoods_father = &likelihoodData_->getD2LikelihoodArray(father->getId());
877  size_t nbSites = _d2Likelihoods_father->size();
878  for (size_t i = 0; i < nbSites; i++)
879  {
880  VVdouble* _d2Likelihoods_father_i = &(*_d2Likelihoods_father)[i];
881  for (size_t c = 0; c < nbClasses_; c++)
882  {
883  Vdouble* _d2Likelihoods_father_i_c = &(*_d2Likelihoods_father_i)[c];
884  for (size_t s = 0; s < nbStates_; s++)
885  {
886  (*_d2Likelihoods_father_i_c)[s] = 1.;
887  }
888  }
889  }
890 
891  size_t nbNodes = father->getNumberOfSons();
892  for (size_t l = 0; l < nbNodes; l++)
893  {
894  const Node* son = father->getSon(l);
895 
896  if (son->getId() == root1_)
897  {
898  const Node* root1 = father->getSon(0);
899  const Node* root2 = father->getSon(1);
900  vector<size_t> * _patternLinks_fatherroot1_ = &likelihoodData_->getArrayPositions(father->getId(), root1->getId());
901  vector<size_t> * _patternLinks_fatherroot2_ = &likelihoodData_->getArrayPositions(father->getId(), root2->getId());
902  VVVdouble* _likelihoodsroot1_ = &likelihoodData_->getLikelihoodArray(root1->getId());
903  VVVdouble* _likelihoodsroot2_ = &likelihoodData_->getLikelihoodArray(root2->getId());
904  double pos = getParameterValue("RootPosition");
905 
906  VVVdouble* d2pxy_root1_ = &d2pxy_[root1_];
907  VVVdouble* d2pxy_root2_ = &d2pxy_[root2_];
908  VVVdouble* dpxy_root1_ = &dpxy_[root1_];
909  VVVdouble* dpxy_root2_ = &dpxy_[root2_];
910  VVVdouble* pxy_root1_ = &pxy_[root1_];
911  VVVdouble* pxy_root2_ = &pxy_[root2_];
912  for (size_t i = 0; i < nbSites; i++)
913  {
914  VVdouble* _likelihoodsroot1__i = &(*_likelihoodsroot1_)[(*_patternLinks_fatherroot1_)[i]];
915  VVdouble* _likelihoodsroot2__i = &(*_likelihoodsroot2_)[(*_patternLinks_fatherroot2_)[i]];
916  VVdouble* _d2Likelihoods_father_i = &(*_d2Likelihoods_father)[i];
917  for (size_t c = 0; c < nbClasses_; c++)
918  {
919  Vdouble* _likelihoodsroot1__i_c = &(*_likelihoodsroot1__i)[c];
920  Vdouble* _likelihoodsroot2__i_c = &(*_likelihoodsroot2__i)[c];
921  Vdouble* _d2Likelihoods_father_i_c = &(*_d2Likelihoods_father_i)[c];
922  VVdouble* d2pxy_root1__c = &(*d2pxy_root1_)[c];
923  VVdouble* d2pxy_root2__c = &(*d2pxy_root2_)[c];
924  VVdouble* dpxy_root1__c = &(*dpxy_root1_)[c];
925  VVdouble* dpxy_root2__c = &(*dpxy_root2_)[c];
926  VVdouble* pxy_root1__c = &(*pxy_root1_)[c];
927  VVdouble* pxy_root2__c = &(*pxy_root2_)[c];
928  for (size_t x = 0; x < nbStates_; x++)
929  {
930  Vdouble* d2pxy_root1__c_x = &(*d2pxy_root1__c)[x];
931  Vdouble* d2pxy_root2__c_x = &(*d2pxy_root2__c)[x];
932  Vdouble* dpxy_root1__c_x = &(*dpxy_root1__c)[x];
933  Vdouble* dpxy_root2__c_x = &(*dpxy_root2__c)[x];
934  Vdouble* pxy_root1__c_x = &(*pxy_root1__c)[x];
935  Vdouble* pxy_root2__c_x = &(*pxy_root2__c)[x];
936  double d2l1 = 0, d2l2 = 0, dl1 = 0, dl2 = 0, l1 = 0, l2 = 0;
937  for (size_t y = 0; y < nbStates_; y++)
938  {
939  d2l1 += (*d2pxy_root1__c_x)[y] * (*_likelihoodsroot1__i_c)[y];
940  d2l2 += (*d2pxy_root2__c_x)[y] * (*_likelihoodsroot2__i_c)[y];
941  dl1 += (*dpxy_root1__c_x)[y] * (*_likelihoodsroot1__i_c)[y];
942  dl2 += (*dpxy_root2__c_x)[y] * (*_likelihoodsroot2__i_c)[y];
943  l1 += (*pxy_root1__c_x)[y] * (*_likelihoodsroot1__i_c)[y];
944  l2 += (*pxy_root2__c_x)[y] * (*_likelihoodsroot2__i_c)[y];
945  }
946  double d2l = pos * pos * d2l1 * l2 + (1. - pos) * (1. - pos) * d2l2 * l1 + 2 * pos * (1. - pos) * dl1 * dl2;
947  (*_d2Likelihoods_father_i_c)[x] *= d2l;
948  }
949  }
950  }
951  }
952  else if (son->getId() == root2_)
953  {
954  //Do nothing, this was accounted when dealing with root1_
955  }
956  else
957  {
958  //Account for a putative multifurcation:
959  vector<size_t> * _patternLinks_father_son = &likelihoodData_->getArrayPositions(father->getId(), son->getId());
960  VVVdouble* _likelihoods_son = &likelihoodData_->getLikelihoodArray(son->getId());
961 
962  VVVdouble* pxy__son = &pxy_[son->getId()];
963  for (size_t i = 0; i < nbSites; i++)
964  {
965  VVdouble* _likelihoods_son_i = &(*_likelihoods_son)[(*_patternLinks_father_son)[i]];
966  VVdouble* _d2Likelihoods_father_i = &(*_d2Likelihoods_father)[i];
967  for (size_t c = 0; c < nbClasses_; c++)
968  {
969  Vdouble* _likelihoods_son_i_c = &(*_likelihoods_son_i)[c];
970  Vdouble* _d2Likelihoods_father_i_c = &(*_d2Likelihoods_father_i)[c];
971  VVdouble* pxy__son_c = &(*pxy__son)[c];
972  for (size_t x = 0; x < nbStates_; x++)
973  {
974  double d2l = 0;
975  Vdouble* pxy__son_c_x = &(*pxy__son_c)[x];
976  for (size_t y = 0; y < nbStates_; y++)
977  {
978  d2l += (*pxy__son_c_x)[y] * (*_likelihoods_son_i_c)[y];
979  }
980  (*_d2Likelihoods_father_i_c)[x] *= d2l;
981  }
982  }
983  }
984  }
985  }
986  return;
987  }
988  else if (variable == "RootPosition")
989  {
990  const Node* father = tree_->getRootNode();
991 
992  // Compute dLikelihoods array for the father node.
993  // Fist initialize to 1:
994  VVVdouble* _d2Likelihoods_father = &likelihoodData_->getD2LikelihoodArray(father->getId());
995  size_t nbSites = _d2Likelihoods_father->size();
996  for (size_t i = 0; i < nbSites; i++)
997  {
998  VVdouble* _d2Likelihoods_father_i = &(*_d2Likelihoods_father)[i];
999  for (size_t c = 0; c < nbClasses_; c++)
1000  {
1001  Vdouble* _d2Likelihoods_father_i_c = &(*_d2Likelihoods_father_i)[c];
1002  for (size_t s = 0; s < nbStates_; s++)
1003  {
1004  (*_d2Likelihoods_father_i_c)[s] = 1.;
1005  }
1006  }
1007  }
1008 
1009  size_t nbNodes = father->getNumberOfSons();
1010  for (size_t l = 0; l < nbNodes; l++)
1011  {
1012  const Node* son = father->getSon(l);
1013 
1014  if (son->getId() == root1_)
1015  {
1016  const Node* root1 = father->getSon(0);
1017  const Node* root2 = father->getSon(1);
1018  vector<size_t> * _patternLinks_fatherroot1_ = &likelihoodData_->getArrayPositions(father->getId(), root1->getId());
1019  vector<size_t> * _patternLinks_fatherroot2_ = &likelihoodData_->getArrayPositions(father->getId(), root2->getId());
1020  VVVdouble* _likelihoodsroot1_ = &likelihoodData_->getLikelihoodArray(root1->getId());
1021  VVVdouble* _likelihoodsroot2_ = &likelihoodData_->getLikelihoodArray(root2->getId());
1022  double len = getParameterValue("BrLenRoot");
1023 
1024  VVVdouble* d2pxy_root1_ = &d2pxy_[root1_];
1025  VVVdouble* d2pxy_root2_ = &d2pxy_[root2_];
1026  VVVdouble* dpxy_root1_ = &dpxy_[root1_];
1027  VVVdouble* dpxy_root2_ = &dpxy_[root2_];
1028  VVVdouble* pxy_root1_ = &pxy_[root1_];
1029  VVVdouble* pxy_root2_ = &pxy_[root2_];
1030  for (size_t i = 0; i < nbSites; i++)
1031  {
1032  VVdouble* _likelihoodsroot1__i = &(*_likelihoodsroot1_)[(*_patternLinks_fatherroot1_)[i]];
1033  VVdouble* _likelihoodsroot2__i = &(*_likelihoodsroot2_)[(*_patternLinks_fatherroot2_)[i]];
1034  VVdouble* _d2Likelihoods_father_i = &(*_d2Likelihoods_father)[i];
1035  for (size_t c = 0; c < nbClasses_; c++)
1036  {
1037  Vdouble* _likelihoodsroot1__i_c = &(*_likelihoodsroot1__i)[c];
1038  Vdouble* _likelihoodsroot2__i_c = &(*_likelihoodsroot2__i)[c];
1039  Vdouble* _d2Likelihoods_father_i_c = &(*_d2Likelihoods_father_i)[c];
1040  VVdouble* d2pxy_root1__c = &(*d2pxy_root1_)[c];
1041  VVdouble* d2pxy_root2__c = &(*d2pxy_root2_)[c];
1042  VVdouble* dpxy_root1__c = &(*dpxy_root1_)[c];
1043  VVdouble* dpxy_root2__c = &(*dpxy_root2_)[c];
1044  VVdouble* pxy_root1__c = &(*pxy_root1_)[c];
1045  VVdouble* pxy_root2__c = &(*pxy_root2_)[c];
1046  for (size_t x = 0; x < nbStates_; x++)
1047  {
1048  Vdouble* d2pxy_root1__c_x = &(*d2pxy_root1__c)[x];
1049  Vdouble* d2pxy_root2__c_x = &(*d2pxy_root2__c)[x];
1050  Vdouble* dpxy_root1__c_x = &(*dpxy_root1__c)[x];
1051  Vdouble* dpxy_root2__c_x = &(*dpxy_root2__c)[x];
1052  Vdouble* pxy_root1__c_x = &(*pxy_root1__c)[x];
1053  Vdouble* pxy_root2__c_x = &(*pxy_root2__c)[x];
1054  double d2l1 = 0, d2l2 = 0, dl1 = 0, dl2 = 0, l1 = 0, l2 = 0;
1055  for (size_t y = 0; y < nbStates_; y++)
1056  {
1057  d2l1 += (*d2pxy_root1__c_x)[y] * (*_likelihoodsroot1__i_c)[y];
1058  d2l2 += (*d2pxy_root2__c_x)[y] * (*_likelihoodsroot2__i_c)[y];
1059  dl1 += (*dpxy_root1__c_x)[y] * (*_likelihoodsroot1__i_c)[y];
1060  dl2 += (*dpxy_root2__c_x)[y] * (*_likelihoodsroot2__i_c)[y];
1061  l1 += (*pxy_root1__c_x)[y] * (*_likelihoodsroot1__i_c)[y];
1062  l2 += (*pxy_root2__c_x)[y] * (*_likelihoodsroot2__i_c)[y];
1063  }
1064  double d2l = len * len * (d2l1 * l2 + d2l2 * l1 - 2 * dl1 * dl2);
1065  (*_d2Likelihoods_father_i_c)[x] *= d2l;
1066  }
1067  }
1068  }
1069  }
1070  else if (son->getId() == root2_)
1071  {
1072  //Do nothing, this was accounted when dealing with root1_
1073  }
1074  else
1075  {
1076  //Account for a putative multifurcation:
1077  vector<size_t> * _patternLinks_father_son = &likelihoodData_->getArrayPositions(father->getId(), son->getId());
1078  VVVdouble* _likelihoods_son = &likelihoodData_->getLikelihoodArray(son->getId());
1079 
1080  VVVdouble* pxy__son = &pxy_[son->getId()];
1081  for (size_t i = 0; i < nbSites; i++)
1082  {
1083  VVdouble* _likelihoods_son_i = &(*_likelihoods_son)[(*_patternLinks_father_son)[i]];
1084  VVdouble* _d2Likelihoods_father_i = &(*_d2Likelihoods_father)[i];
1085  for (size_t c = 0; c < nbClasses_; c++)
1086  {
1087  Vdouble* _likelihoods_son_i_c = &(*_likelihoods_son_i)[c];
1088  Vdouble* _d2Likelihoods_father_i_c = &(*_d2Likelihoods_father_i)[c];
1089  VVdouble* pxy__son_c = &(*pxy__son)[c];
1090  for (size_t x = 0; x < nbStates_; x++)
1091  {
1092  double d2l = 0;
1093  Vdouble* pxy__son_c_x = &(*pxy__son_c)[x];
1094  for (size_t y = 0; y < nbStates_; y++)
1095  {
1096  d2l += (*pxy__son_c_x)[y] * (*_likelihoods_son_i_c)[y];
1097  }
1098  (*_d2Likelihoods_father_i_c)[x] *= d2l;
1099  }
1100  }
1101  }
1102  }
1103  }
1104  return;
1105  }
1106 
1107  // Get the node with the branch whose length must be derivated:
1108  size_t brI = TextTools::to<size_t>(variable.substr(5));
1109  const Node* branch = nodes_[brI];
1110  const Node* father = branch->getFather();
1111 
1112  // Compute dLikelihoods array for the father node.
1113  // Fist initialize to 1:
1114  VVVdouble* _d2Likelihoods_father = &likelihoodData_->getD2LikelihoodArray(father->getId());
1115  size_t nbSites = _d2Likelihoods_father->size();
1116  for (size_t i = 0; i < nbSites; i++)
1117  {
1118  VVdouble* _d2Likelihoods_father_i = &(*_d2Likelihoods_father)[i];
1119  for (size_t c = 0; c < nbClasses_; c++)
1120  {
1121  Vdouble* _d2Likelihoods_father_i_c = &(*_d2Likelihoods_father_i)[c];
1122  for (size_t s = 0; s < nbStates_; s++)
1123  {
1124  (*_d2Likelihoods_father_i_c)[s] = 1.;
1125  }
1126  }
1127  }
1128 
1129  size_t nbNodes = father->getNumberOfSons();
1130  for (size_t l = 0; l < nbNodes; l++)
1131  {
1132  const Node* son = father->getSon(l);
1133 
1134  vector<size_t> * _patternLinks_father_son = &likelihoodData_->getArrayPositions(father->getId(), son->getId());
1135  VVVdouble* _likelihoods_son = &likelihoodData_->getLikelihoodArray(son->getId());
1136 
1137  if (son == branch)
1138  {
1139  VVVdouble* d2pxy__son = &d2pxy_[son->getId()];
1140  for (size_t i = 0; i < nbSites; i++)
1141  {
1142  VVdouble* _likelihoods_son_i = &(*_likelihoods_son)[(*_patternLinks_father_son)[i]];
1143  VVdouble* _d2Likelihoods_father_i = &(*_d2Likelihoods_father)[i];
1144  for (size_t c = 0; c < nbClasses_; c++)
1145  {
1146  Vdouble* _likelihoods_son_i_c = &(*_likelihoods_son_i)[c];
1147  Vdouble* _d2Likelihoods_father_i_c = &(*_d2Likelihoods_father_i)[c];
1148  VVdouble* d2pxy__son_c = &(*d2pxy__son)[c];
1149  for (size_t x = 0; x < nbStates_; x++)
1150  {
1151  double d2l = 0;
1152  Vdouble* d2pxy__son_c_x = &(*d2pxy__son_c)[x];
1153  for (size_t y = 0; y < nbStates_; y++)
1154  {
1155  d2l += (*d2pxy__son_c_x)[y] * (*_likelihoods_son_i_c)[y];
1156  }
1157  (*_d2Likelihoods_father_i_c)[x] *= d2l;
1158  }
1159  }
1160  }
1161  }
1162  else
1163  {
1164  VVVdouble* pxy__son = &pxy_[son->getId()];
1165  for (size_t i = 0; i < nbSites; i++)
1166  {
1167  VVdouble* _likelihoods_son_i = &(*_likelihoods_son)[(*_patternLinks_father_son)[i]];
1168  VVdouble* _d2Likelihoods_father_i = &(*_d2Likelihoods_father)[i];
1169  for (size_t c = 0; c < nbClasses_; c++)
1170  {
1171  Vdouble* _likelihoods_son_i_c = &(*_likelihoods_son_i)[c];
1172  Vdouble* _d2Likelihoods_father_i_c = &(*_d2Likelihoods_father_i)[c];
1173  VVdouble* pxy__son_c = &(*pxy__son)[c];
1174  for (size_t x = 0; x < nbStates_; x++)
1175  {
1176  double d2l = 0;
1177  Vdouble* pxy__son_c_x = &(*pxy__son_c)[x];
1178  for (size_t y = 0; y < nbStates_; y++)
1179  {
1180  d2l += (*pxy__son_c_x)[y] * (*_likelihoods_son_i_c)[y];
1181  }
1182  (*_d2Likelihoods_father_i_c)[x] *= d2l;
1183  }
1184  }
1185  }
1186  }
1187  }
1188 
1189  // Now we go down the tree toward the root node:
1191 }
1192 
1193 /******************************************************************************/
1194 
1196 {
1197  const Node* father = node->getFather();
1198  // We assume that the _dLikelihoods array has been filled for the current node 'node'.
1199  // We will evaluate the array for the father node.
1200  if (father == 0) return; // We reached the root!
1201 
1202  // Compute dLikelihoods array for the father node.
1203  // Fist initialize to 1:
1204  VVVdouble* _d2Likelihoods_father = &likelihoodData_->getD2LikelihoodArray(father->getId());
1205  size_t nbSites = _d2Likelihoods_father->size();
1206  for (size_t i = 0; i < nbSites; i++)
1207  {
1208  VVdouble* _d2Likelihoods_father_i = &(*_d2Likelihoods_father)[i];
1209  for (size_t c = 0; c < nbClasses_; c++)
1210  {
1211  Vdouble* _d2Likelihoods_father_i_c = &(*_d2Likelihoods_father_i)[c];
1212  for (size_t s = 0; s < nbStates_; s++)
1213  {
1214  (*_d2Likelihoods_father_i_c)[s] = 1.;
1215  }
1216  }
1217  }
1218 
1219  size_t nbNodes = father->getNumberOfSons();
1220  for (size_t l = 0; l < nbNodes; l++)
1221  {
1222  const Node* son = father->getSon(l);
1223 
1224  VVVdouble* pxy__son = &pxy_[son->getId()];
1225  vector<size_t> * _patternLinks_father_son = &likelihoodData_->getArrayPositions(father->getId(), son->getId());
1226 
1227  if (son == node)
1228  {
1229  VVVdouble* _d2Likelihoods_son = &likelihoodData_->getD2LikelihoodArray(son->getId());
1230  for (size_t i = 0; i < nbSites; i++)
1231  {
1232  VVdouble* _d2Likelihoods_son_i = &(*_d2Likelihoods_son)[(*_patternLinks_father_son)[i]];
1233  VVdouble* _d2Likelihoods_father_i = &(*_d2Likelihoods_father)[i];
1234  for (size_t c = 0; c < nbClasses_; c++)
1235  {
1236  Vdouble* _d2Likelihoods_son_i_c = &(*_d2Likelihoods_son_i)[c];
1237  Vdouble* _d2Likelihoods_father_i_c = &(*_d2Likelihoods_father_i)[c];
1238  VVdouble* pxy__son_c = &(*pxy__son)[c];
1239  for (size_t x = 0; x < nbStates_; x++)
1240  {
1241  double d2l = 0;
1242  Vdouble* pxy__son_c_x = &(*pxy__son_c)[x];
1243  for (size_t y = 0; y < nbStates_; y++)
1244  {
1245  d2l += (*pxy__son_c_x)[y] * (*_d2Likelihoods_son_i_c)[y];
1246  }
1247  (*_d2Likelihoods_father_i_c)[x] *= d2l;
1248  }
1249  }
1250  }
1251  }
1252  else
1253  {
1254  VVVdouble* _likelihoods_son = &likelihoodData_->getLikelihoodArray(son->getId());
1255  for (size_t i = 0; i < nbSites; i++)
1256  {
1257  VVdouble* _likelihoods_son_i = &(*_likelihoods_son)[(*_patternLinks_father_son)[i]];
1258  VVdouble* _d2Likelihoods_father_i = &(*_d2Likelihoods_father)[i];
1259  for (size_t c = 0; c < nbClasses_; c++)
1260  {
1261  Vdouble* _likelihoods_son_i_c = &(*_likelihoods_son_i)[c];
1262  Vdouble* _d2Likelihoods_father_i_c = &(*_d2Likelihoods_father_i)[c];
1263  VVdouble* pxy__son_c = &(*pxy__son)[c];
1264  for (size_t x = 0; x < nbStates_; x++)
1265  {
1266  double dl = 0;
1267  Vdouble* pxy__son_c_x = &(*pxy__son_c)[x];
1268  for (size_t y = 0; y < nbStates_; y++)
1269  {
1270  dl += (*pxy__son_c_x)[y] * (*_likelihoods_son_i_c)[y];
1271  }
1272  (*_d2Likelihoods_father_i_c)[x] *= dl;
1273  }
1274  }
1275  }
1276  }
1277  }
1278 
1279  //Next step: move toward grand father...
1281 }
1282 
1283 /******************************************************************************/
1284 
1286 {
1287  computeSubtreeLikelihood(tree_->getRootNode());
1288 }
1289 
1290 /******************************************************************************/
1291 
1293 {
1294  if (node->isLeaf()) return;
1295 
1296  size_t nbSites = likelihoodData_->getLikelihoodArray(node->getId()).size();
1297  size_t nbNodes = node->getNumberOfSons();
1298 
1299  // Must reset the likelihood array first (i.e. set all of them to 1):
1300  VVVdouble* _likelihoods_node = &likelihoodData_->getLikelihoodArray(node->getId());
1301  for (size_t i = 0; i < nbSites; i++)
1302  {
1303  //For each site in the sequence,
1304  VVdouble* _likelihoods_node_i = &(*_likelihoods_node)[i];
1305  for (size_t c = 0; c < nbClasses_; c++)
1306  {
1307  //For each rate classe,
1308  Vdouble* _likelihoods_node_i_c = &(*_likelihoods_node_i)[c];
1309  for (size_t x = 0; x < nbStates_; x++)
1310  {
1311  //For each initial state,
1312  (*_likelihoods_node_i_c)[x] = 1.;
1313  }
1314  }
1315  }
1316 
1317  for (size_t l = 0; l < nbNodes; l++)
1318  {
1319  //For each son node,
1320 
1321  const Node* son = node->getSon(l);
1322 
1323  computeSubtreeLikelihood(son); //Recursive method:
1324 
1325  VVVdouble* pxy__son = &pxy_[son->getId()];
1326  vector<size_t> * _patternLinks_node_son = &likelihoodData_->getArrayPositions(node->getId(), son->getId());
1327  VVVdouble* _likelihoods_son = &likelihoodData_->getLikelihoodArray(son->getId());
1328 
1329  for (size_t i = 0; i < nbSites; i++)
1330  {
1331  //For each site in the sequence,
1332  VVdouble* _likelihoods_son_i = &(*_likelihoods_son)[(*_patternLinks_node_son)[i]];
1333  VVdouble* _likelihoods_node_i = &(*_likelihoods_node)[i];
1334  for (size_t c = 0; c < nbClasses_; c++)
1335  {
1336  //For each rate classe,
1337  Vdouble* _likelihoods_son_i_c = &(*_likelihoods_son_i)[c];
1338  Vdouble* _likelihoods_node_i_c = &(*_likelihoods_node_i)[c];
1339  VVdouble* pxy__son_c = &(*pxy__son)[c];
1340  for (size_t x = 0; x < nbStates_; x++)
1341  {
1342  //For each initial state,
1343  Vdouble* pxy__son_c_x = &(*pxy__son_c)[x];
1344  double likelihood = 0;
1345  for (size_t y = 0; y < nbStates_; y++)
1346  {
1347  likelihood += (*pxy__son_c_x)[y] * (*_likelihoods_son_i_c)[y];
1348  }
1349  (*_likelihoods_node_i_c)[x] *= likelihood;
1350  }
1351  }
1352  }
1353  }
1354 }
1355 
1356 
1357 /******************************************************************************/
1358 
1360 {
1361  cout << "Likelihoods at node " << node->getName() << ": " << endl;
1363  cout << " ***" << endl;
1364 }
1365 
virtual void computeDownSubtreeDLikelihood(const Node *)
This class implement the &#39;traditional&#39; way of computing likelihood for a tree, allowing for non-homog...
Substitution models manager for non-homogeneous / non-reversible models of evolution.
Partial implementation for branch non-homogeneous models of the TreeLikelihood interface.
virtual void computeTreeD2Likelihood(const std::string &variable)
std::vector< Node * > nodes_
Pointer toward all nodes in the tree.
virtual std::string getName() const
Get the name associated to this node, if there is one, otherwise throw a NodeException.
Definition: Node.h:236
size_t getRootArrayPosition(size_t currentPosition) const
void setTree(const TreeTemplate< Node > *tree)
Set the tree associated to the data.
virtual void computeAllTransitionProbabilities()
Fill the pxy_, dpxy_ and d2pxy_ arrays for all nodes.
double getLikelihoodForASite(size_t site) const
Get the likelihood for a site.
virtual void applyParameters()
All parameters are stored in a parameter list. This function apply these parameters to the substituti...
virtual const Node * getSon(size_t pos) const
Definition: Node.h:395
std::vector< double > getRootFrequencies() const
STL namespace.
double getFirstOrderDerivative(const std::string &variable) const
void fireParameterChanged(const ParameterList &params)
discrete Rate Across Sites, (simple) Recursive likelihood data structure.
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
RNonHomogeneousTreeLikelihood & operator=(const RNonHomogeneousTreeLikelihood &lik)
static void displayLikelihoodArray(const VVVdouble &likelihoodArray)
Print the likelihood array to terminal (debugging tool).
virtual bool isLeaf() const
Definition: Node.h:692
ParameterList getNodeParameters() const
Get the parameters corresponding attached to the nodes of the tree.
AbstractNonHomogeneousTreeLikelihood & operator=(const AbstractNonHomogeneousTreeLikelihood &lik)
Assignation operator.
VVVdouble & getDLikelihoodArray(int nodeId)
double getLogLikelihoodForASite(size_t site) const
Get the logarithm of the likelihood for a site.
virtual int getId() const
Get the node&#39;s id.
Definition: Node.h:203
virtual double getD2LikelihoodForASiteForARateClass(size_t site, size_t rateClass) const
double getLikelihoodForASiteForARateClass(size_t site, size_t rateClass) const
Get the likelihood for a site knowing its rate class.
virtual double getD2LikelihoodForASite(size_t site) const
void setData(const SiteContainer &sites)
Set the dataset for which the likelihood must be evaluated.
static SiteContainer * getSequenceSubset(const SiteContainer &sequenceSet, const Node &node)
Extract the sequences corresponding to a given subtree.
RNonHomogeneousTreeLikelihood(const Tree &tree, SubstitutionModelSet *modelSet, DiscreteDistribution *rDist, bool verbose=true, bool usePatterns=true, bool reparametrizeRoot=false)
Build a new NonHomogeneousTreeLikelihood object without data.
virtual void computeTreeDLikelihood(const std::string &variable)
virtual void computeTransitionProbabilitiesForNode(const Node *node)
Fill the pxy_, dpxy_ and d2pxy_ arrays for one node.
virtual void computeSubtreeLikelihood(const Node *node)
Compute the likelihood for a subtree defined by the Tree::Node node.
virtual double getDLikelihoodForASite(size_t site) const
std::map< int, const Node * > idToNode_
An index linking nodes to their id, for faster access than the getNode() method.
virtual double getDLikelihoodForASiteForARateClass(size_t site, size_t rateClass) const
virtual void displayLikelihood(const Node *node)
This method is mainly for debugging purpose.
The phylogenetic node class.
Definition: Node.h:90
virtual double getDLogLikelihoodForASite(size_t site) const
double getSecondOrderDerivative(const std::string &variable) const
void init_(bool usePatterns)
Method called by constructors.
DRASRTreeLikelihoodData * clone() const
void setParameters(const ParameterList &parameters)
Implements the Function interface.
double getLikelihood() const
Get the likelihood for the whole dataset.
double getLogLikelihood() const
Get the logarithm of the likelihood for the whole dataset.
virtual size_t getNumberOfSons() const
Definition: Node.h:388
virtual double getD2LogLikelihoodForASite(size_t site) const
const std::vector< size_t > & getArrayPositions(int parentId, int sonId) const
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 getLogLikelihoodForASiteForARateClass(size_t site, size_t rateClass) const
Get the logarithm of the likelihood for a site knowing its rate class.
VVVdouble & getLikelihoodArray(int nodeId)
std::vector< int > getNodesWithParameter(const std::string &name) const
VVVdouble & getD2LikelihoodArray(int nodeId)
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...