bpp-phyl  2.1.0
 All Classes Namespaces Files Functions Variables Friends Pages
AbstractSubstitutionModel.cpp
Go to the documentation of this file.
1 //
2 // File: AbstractSubstitutionModel.cpp
3 // Created by: Julien Dutheil
4 // Created on: Tue May 27 10:31:49 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 
42 #include <Bpp/Text/TextTools.h>
47 
48 // From SeqLib:
50 
51 using namespace bpp;
52 using namespace std;
53 
54 /******************************************************************************/
55 
56 AbstractSubstitutionModel::AbstractSubstitutionModel(const Alphabet* alpha, const std::string& prefix) :
58  alphabet_(alpha),
59  size_(alpha->getSize()),
60  rate_(1),
61  chars_(size_),
62  generator_(size_, size_),
63  freq_(size_),
64  exchangeability_(size_, size_),
65  pijt_(size_, size_),
66  dpijt_(size_, size_),
67  d2pijt_(size_, size_),
68  eigenDecompose_(true),
69  eigenValues_(size_),
70  iEigenValues_(size_),
71  isDiagonalizable_(false),
72  rightEigenVectors_(size_, size_),
73  isNonSingular_(false),
74  leftEigenVectors_(size_, size_),
75  vPowGen_(),
76  tmpMat_(size_, size_)
77 {
78  for (size_t i = 0; i < size_; i++)
79  {
80  freq_[i] = 1.0 / static_cast<double>(size_);
81  chars_[i] = static_cast<int>(i);
82  }
83 }
84 
85 /******************************************************************************/
86 
88 {
89  // if the object is not an AbstractReversibleSubstitutionModel,
90  // computes the exchangeability_ Matrix (otherwise the generator_
91  // has been computed from the exchangeability_)
92 
93  if (!dynamic_cast<AbstractReversibleSubstitutionModel*>(this)) {
94  for (size_t i = 0; i < size_; i++)
95  {
96  for (size_t j = 0; j < size_; j++)
97  {
98  exchangeability_(i, j) = generator_(i, j) / freq_[j];
99  }
100  }
101  }
102 
103  // Compute eigen values and vectors:
105  {
107  rightEigenVectors_ = ev.getV();
110  try
111  {
113  isNonSingular_ = true;
114  isDiagonalizable_ = true;
115  for (size_t i = 0; i < size_ && isDiagonalizable_; i++)
116  {
117  if (abs(iEigenValues_[i]) > NumConstants::TINY())
118  isDiagonalizable_ = false;
119  }
120  }
121  catch (ZeroDivisionException& e)
122  {
123  ApplicationTools::displayMessage("Singularity during diagonalization. Taylor series used instead.");
124 
125  isNonSingular_ = false;
126  isDiagonalizable_ = false;
128  }
129  }
130 }
131 
132 
133 /******************************************************************************/
134 
136 {
137  if (t == 0)
138  {
140  }
141  else if (isNonSingular_)
142  {
143  if (isDiagonalizable_)
144  {
145  MatrixTools::mult<double>(rightEigenVectors_, VectorTools::exp(eigenValues_ * (rate_ * t)), leftEigenVectors_, pijt_);
146  }
147  else
148  {
149  std::vector<double> vdia(size_);
150  std::vector<double> vup(size_ - 1);
151  std::vector<double> vlo(size_ - 1);
152  double c = 0, s = 0;
153  double l = rate_ * t;
154  for (size_t i = 0; i < size_; i++)
155  {
156  vdia[i] = std::exp(eigenValues_[i] * l);
157  if (iEigenValues_[i] != 0)
158  {
159  s = std::sin(iEigenValues_[i] * l);
160  c = std::cos(iEigenValues_[i] * l);
161  vdia[i] *= c;
162  vup[i] = vdia[i] * s;
163  vlo[i] = -vup[i];
164  vdia[i + 1] = vdia[i]; // trick to avoid computation
165  i++;
166  }
167  else
168  {
169  if (i < size_ - 1)
170  {
171  vup[i] = 0;
172  vlo[i] = 0;
173  }
174  }
175  }
176  MatrixTools::mult<double>(rightEigenVectors_, vdia, vup, vlo, leftEigenVectors_, pijt_);
177  }
178  }
179  else
180  {
182  double s = 1.0;
183  double v = rate_ * t;
184  size_t m = 0;
185  while (v > 0.5) // exp(r*t*A)=(exp(r*t/(2^m) A))^(2^m)
186  {
187  m += 1;
188  v /= 2;
189  }
190  for (size_t i = 1; i < vPowGen_.size(); i++)
191  {
192  s *= v / static_cast<double>(i);
194  }
195  while (m > 0) // recover the 2^m
196  {
199  m--;
200  }
201  }
202  //MatrixTools::print(pijt_);
203  return pijt_;
204 }
205 
207 {
208  if (isNonSingular_)
209  {
210  if (isDiagonalizable_)
211  {
213  }
214  else
215  {
216  std::vector<double> vdia(size_);
217  std::vector<double> vup(size_ - 1);
218  std::vector<double> vlo(size_ - 1);
219  double c, s, e;
220  double l = rate_ * t;
221  for (size_t i = 0; i < size_; i++)
222  {
223  e = std::exp(eigenValues_[i] * l);
224  if (iEigenValues_[i] != 0)
225  {
226  s = std::sin(iEigenValues_[i] * l);
227  c = std::cos(iEigenValues_[i] * l);
228  vdia[i] = rate_ * (eigenValues_[i] * c - iEigenValues_[i] * s) * e;
229  vup[i] = rate_ * (eigenValues_[i] * s + iEigenValues_[i] * c) * e;
230  vlo[i] = -vup[i];
231  vdia[i + 1] = vdia[i]; // trick to avoid computation
232  i++;
233  }
234  else
235  {
236  if (i < size_ - 1)
237  {
238  vup[i] = 0;
239  vlo[i] = 0;
240  }
241  }
242  }
243  MatrixTools::mult<double>(rightEigenVectors_, vdia, vup, vlo, leftEigenVectors_, dpijt_);
244  }
245  }
246  else
247  {
249  double s = 1.0;
250  double v = rate_ * t;
251  size_t m = 0;
252  while (v > 0.5) // r*A*exp(t*r*A)=r*A*(exp(r*t/(2^m) A))^(2^m)
253  {
254  m += 1;
255  v /= 2;
256  }
257  for (size_t i = 1; i < vPowGen_.size(); i++)
258  {
259  s *= v / static_cast<double>(i);
261  }
262  while (m > 0) // recover the 2^m
263  {
266  m--;
267  }
271  }
272  return dpijt_;
273 }
274 
276 {
277  if (isNonSingular_)
278  {
279  if (isDiagonalizable_)
280  {
282  }
283  else
284  {
285  std::vector<double> vdia(size_);
286  std::vector<double> vup(size_ - 1);
287  std::vector<double> vlo(size_ - 1);
288  double c, s, e;
289  double l = rate_ * t;
290  for (size_t i = 0; i < size_; i++)
291  {
292  e = std::exp(eigenValues_[i] * l);
293  if (iEigenValues_[i] != 0)
294  {
295  s = std::sin(iEigenValues_[i] * l);
296  c = std::cos(iEigenValues_[i] * l);
297  vdia[i] = NumTools::sqr(rate_)
299  - 2 * eigenValues_[i] * iEigenValues_[i] * s) * e;
300  vup[i] = NumTools::sqr(rate_)
302  - 2 * eigenValues_[i] * iEigenValues_[i] * c) * e;
303  vlo[i] = -vup[i];
304  vdia[i + 1] = vdia[i]; // trick to avoid computation
305  i++;
306  }
307  else
308  {
309  if (i < size_ - 1)
310  {
311  vup[i] = 0;
312  vlo[i] = 0;
313  }
314  }
315  }
316  MatrixTools::mult<double>(rightEigenVectors_, vdia, vup, vlo, leftEigenVectors_, d2pijt_);
317  }
318  }
319  else
320  {
322  double s = 1.0;
323  double v = rate_ * t;
324  size_t m = 0;
325  while (v > 0.5) // r^2*A^2*exp(t*r*A)=r^2*A^2*(exp(r*t/(2^m) A))^(2^m)
326  {
327  m += 1;
328  v /= 2;
329  }
330  for (size_t i = 1; i < vPowGen_.size(); i++)
331  {
332  s *= v / static_cast<double>(i);
334  }
335  while (m > 0) // recover the 2^m
336  {
339  m--;
340  }
344  }
345  return d2pijt_;
346 }
347 
348 /******************************************************************************/
349 
351 {
352  if (i >= size_)
353  throw IndexOutOfBoundsException("AbstractSubstitutionModel::getInitValue", i, 0, size_ - 1);
354  if (state < 0 || !alphabet_->isIntInAlphabet(state))
355  throw BadIntException(state, "AbstractSubstitutionModel::getInitValue. Character " + alphabet_->intToChar(state) + " is not allowed in model.");
356  vector<int> states = alphabet_->getAlias(state);
357  for (size_t j = 0; j < states.size(); j++)
358  {
359  if (getAlphabetChar(i) == states[j])
360  return 1.;
361  }
362  return 0.;
363 }
364 
365 /******************************************************************************/
366 
368 {
369  map<int, int> counts;
370  SequenceContainerTools::getCounts(data, counts);
371  double t = 0;
372  map<int, double> freqs;
373 
374  for (int i = 0; i < static_cast<int>(size_); i++)
375  {
376  t += counts[i] + pseudoCount;
377  }
378  for (int i = 0; i < static_cast<int>(size_); i++)
379  {
380  freqs[i] = (static_cast<double>(counts[i]) + pseudoCount) / t;
381  }
382 
383  // Re-compute generator and eigen values:
384  setFreq(freqs);
385 }
386 
387 /******************************************************************************/
388 
389 void AbstractSubstitutionModel::setFreq(map<int, double>& freqs)
390 {
391  for (int i = 0; i < static_cast<int>(size_); i++)
392  {
393  freq_[i] = freqs[i];
394  }
395  // Re-compute generator and eigen values:
396  updateMatrices();
397 }
398 
399 /******************************************************************************/
400 
402 {
403  vector<double> v;
405  return -VectorTools::scalar<double, double>(v, freq_);
406 }
407 
408 /******************************************************************************/
409 
412 }
413 
414 /******************************************************************************/
415 
417 {
418  return rate_;
419 }
420 
421 /******************************************************************************/
422 
424 {
425  if (rate <= 0)
426  throw Exception("Bad value for rate: " + TextTools::toString(rate));
427 
428  if (hasParameter("rate"))
429  setParameterValue("rate", rate_);
430 
431  rate_ = rate;
432 }
433 
435 {
436  addParameter_(new Parameter(getNamespace() + "rate", rate_, &Parameter::R_PLUS_STAR));
437 }
438 
439 /******************************************************************************/
440 
442 {
445  MatrixTools::mult(exchangeability_, Pi, generator_); // Diagonal elements of the exchangability matrix will be ignored.
446  // Compute diagonal elements of the generator:
447  for (size_t i = 0; i < size_; i++)
448  {
449  double lambda = 0;
450  for (size_t j = 0; j < size_; j++)
451  {
452  if (j != i)
453  lambda += generator_(i, j);
454  }
455  generator_(i, i) = -lambda;
456  }
457  // Normalization:
458  double scale = getScale();
459  MatrixTools::scale(generator_, 1. / scale);
460 
461  // Normalize exchangeability matrix too:
463  // Compute diagonal elements of the exchangeability matrix:
464  for (size_t i = 0; i < size_; i++)
465  {
466  exchangeability_(i, i) = generator_(i, i) / freq_[i];
467  }
469 }
470 
471 /******************************************************************************/
472