bpp-phyl3 3.0.0
YpR.cpp
Go to the documentation of this file.
1// SPDX-FileCopyrightText: The Bio++ Development Group
2//
3// SPDX-License-Identifier: CECILL-2.1
4
5#include "YpR.h"
6
7// From the STL:
8#include <cmath>
9
10using namespace bpp;
11
15
17
18#include <Bpp/Text/TextTools.h>
19
20/******************************************************************************/
21
22YpR::YpR(
23 std::shared_ptr<const RNY> alph,
24 std::unique_ptr<NucleotideSubstitutionModelInterface> pm,
25 const string& prefix) :
27 AbstractSubstitutionModel(alph, make_shared<CanonicalStateMap>(alph, false), prefix),
28 pmodel_(pm->clone()),
29 nestedPrefix_(pm->getNamespace())
30{
31 pmodel_->setNamespace(prefix + nestedPrefix_);
32 pmodel_->enableEigenDecomposition(0);
33 pmodel_->computeFrequencies(false);
34
35 addParameters_(pmodel_->getParameters());
37}
38
39YpR::YpR(const YpR& ypr) :
42 pmodel_(ypr.pmodel_->clone()),
43 nestedPrefix_(ypr.getNestedPrefix())
44{}
45
47{
48 updateMatrices_(0, 0, 0, 0, 0, 0, 0, 0);
49}
50
51
52void YpR::updateMatrices_(double CgT, double cGA,
53 double TgC, double tGA,
54 double CaT, double cAG,
55 double TaC, double tAC)
56{
58
59 // Generator:
60 unsigned int i, j, i1, i2, i3, j1, j2, j3;
61
62 std::vector<double> a(4); // a[A], a[G], a[C], a[T]
63 std::vector<double> b(4); // b[A], b[G], b[C], b[T]
64
65 // From M
66
67 a[0] = pmodel_->Qij(2, 0);
68 a[1] = pmodel_->Qij(0, 2);
69 a[2] = pmodel_->Qij(3, 1);
70 a[3] = pmodel_->Qij(1, 3);
71 b[0] = (pmodel_->Qij(1, 0) + pmodel_->Qij(3, 0)) / 2; // To limit numerical issues
72 b[1] = (pmodel_->Qij(1, 2) + pmodel_->Qij(3, 2)) / 2;
73 b[2] = (pmodel_->Qij(0, 1) + pmodel_->Qij(2, 1)) / 2;
74 b[3] = (pmodel_->Qij(0, 3) + pmodel_->Qij(2, 3)) / 2;
75
76
77 // M_1 on R C T
78 RowMatrix<double> M1(3, 3);
79
80 M1(0, 0) = 0;
81 M1(0, 1) = b[2];
82 M1(0, 2) = b[3];
83 M1(1, 0) = pmodel_->Qij(1, 0) + pmodel_->Qij(1, 2);
84 M1(1, 1) = 0;
85 M1(1, 2) = pmodel_->Qij(1, 3);
86 M1(2, 0) = pmodel_->Qij(3, 0) + pmodel_->Qij(3, 2);
87 M1(2, 1) = pmodel_->Qij(3, 1);
88 M1(2, 2) = 0;
89
90 // M_2 on A G C T
91 RowMatrix<double> M2(4, 4);
92
93 M2(0, 0) = 0;
94 M2(0, 1) = pmodel_->Qij(0, 2);
95 M2(0, 2) = pmodel_->Qij(0, 1);
96 M2(0, 3) = pmodel_->Qij(0, 3);
97 M2(1, 0) = pmodel_->Qij(2, 0);
98 M2(1, 1) = 0;
99 M2(1, 2) = pmodel_->Qij(2, 1);
100 M2(1, 3) = pmodel_->Qij(2, 3);
101 M2(2, 0) = pmodel_->Qij(1, 0);
102 M2(2, 1) = pmodel_->Qij(1, 2);
103 M2(2, 2) = 0;
104 M2(2, 3) = pmodel_->Qij(1, 3);
105 M2(3, 0) = pmodel_->Qij(3, 0);
106 M2(3, 1) = pmodel_->Qij(3, 2);
107 M2(3, 2) = pmodel_->Qij(3, 1);
108 M2(3, 3) = 0;
109
110 // M_3 on A G Y
111 RowMatrix<double> M3(3, 3);
112
113 M3(0, 0) = 0;
114 M3(0, 1) = pmodel_->Qij(0, 2);
115 M3(0, 2) = pmodel_->Qij(0, 1) + pmodel_->Qij(0, 3);
116 M3(1, 0) = pmodel_->Qij(2, 0);
117 M3(1, 1) = 0;
118 M3(1, 2) = pmodel_->Qij(2, 1) + pmodel_->Qij(2, 3);
119 M3(2, 0) = b[0];
120 M3(2, 1) = b[1];
121 M3(2, 2) = 0;
122
123
124 for (i1 = 0; i1 < 3; i1++)
125 {
126 for (i2 = 0; i2 < 4; i2++)
127 {
128 for (i3 = 0; i3 < 3; i3++)
129 {
130 i = 12 * i1 + 3 * i2 + i3;
131 for (j1 = 0; j1 < 3; j1++)
132 {
133 for (j2 = 0; j2 < 4; j2++)
134 {
135 for (j3 = 0; j3 < 3; j3++)
136 {
137 j = 12 * j1 + 3 * j2 + j3;
138 if ((i1 == j1) && (i2 == j2))
139 generator_(i, j) = M3(i3, j3);
140 else if ((i1 == j1) && (i3 == j3))
141 generator_(i, j) = M2(i2, j2);
142 else if ((i2 == j2) && (i3 == j3))
143 generator_(i, j) = M1(i1, j1);
144 else
145 generator_(i, j) = 0;
146 }
147 }
148 }
149 }
150 }
151 }
152
153 // Introduction des dependances
154
155 for (i3 = 0; i3 < 3; i3++)
156 {
157 generator_(15 + i3, 12 + i3) += cGA * a[0]; // CG -> CA
158 generator_(12 * i3 + 7, 12 * i3 + 6) += cGA * a[0];
159
160 generator_(15 + i3, 27 + i3) += CgT * a[3]; // CG -> TG
161 generator_(12 * i3 + 7, 12 * i3 + 10) += CgT * a[3];
162
163 generator_(27 + i3, 24 + i3) += tGA * a[0]; // TG -> TA
164 generator_(12 * i3 + 10, 12 * i3 + 9) += tGA * a[0];
165
166 generator_(27 + i3, 15 + i3) += TgC * a[2]; // TG -> CG
167 generator_(12 * i3 + 10, 12 * i3 + 7) += TgC * a[2];
168
169 generator_(12 + i3, 24 + i3) += CaT * a[3]; // CA -> TA
170 generator_(12 * i3 + 6, 12 * i3 + 9) += CaT * a[3];
171
172 generator_(12 + i3, 15 + i3) += cAG * a[1]; // CA -> CG
173 generator_(12 * i3 + 6, 12 * i3 + 7) += cAG * a[1];
174
175 generator_(24 + i3, 27 + i3) += tAC * a[1]; // TA -> TG
176 generator_(12 * i3 + 9, 12 * i3 + 10) += tAC * a[1];
177
178 generator_(24 + i3, 12 + i3) += TaC * a[2]; // TA -> CA
179 generator_(12 * i3 + 9, 12 * i3 + 6) += TaC * a[2];
180 }
181
182 double x;
183
184 for (i = 0; i < 36; ++i)
185 {
186 x = 0;
187 for (j = 0; j < 36; ++j)
188 {
189 if (j != i)
190 x += generator_(i, j);
191 }
192 generator_(i, i) = -x;
193 }
194
195 // calcul spectral
196
200
202
203 try
204 {
206
207 isNonSingular_ = true;
208 isDiagonalizable_ = true;
209 for (i = 0; i < size_; i++)
210 {
212 {
213 isDiagonalizable_ = false;
214 }
215 }
216
217 // frequence stationnaire
218
219 x = 0;
220 j = 0;
221 while (j < 36)
222 {
225 {
226 eigenValues_[j] = 0; // to avoid approximation problems in the future
227 for (i = 0; i < 36; i++)
228 {
229 freq_[i] = leftEigenVectors_(j, i);
230 x += freq_[i];
231 }
232 break;
233 }
234 j++;
235 }
236 for (i = 0; i < 36; i++)
237 {
238 freq_[i] /= x;
239 }
240 }
241 catch (ZeroDivisionException& e)
242 {
243 ApplicationTools::displayMessage("Singularity during diagonalization. Taylor series used instead.");
244 isNonSingular_ = false;
245 isDiagonalizable_ = false;
246
247 if (vPowGen_.size() == 0)
248 vPowGen_.resize(30);
249
250 double min = generator_(0, 0);
251 for (i = 1; i < 36; i++)
252 {
253 if (min > generator_(i, i))
254 min = generator_(i, i);
255 }
256
257 setScale(-1 / min);
258
259 MatrixTools::getId(36, tmpMat_); // to compute the equilibrium frequency (Q+Id)^256
260
263
264 for (i = 0; i < 36; i++)
265 {
266 freq_[i] = vPowGen_[0](0, i);
267 }
268
270 }
271
272 // mise a l'echelle
273
274 x = 0;
275 for (i1 = 0; i1 < 3; i1++)
276 {
277 for (i2 = 0; i2 < 4; i2++)
278 {
279 for (i3 = 0; i3 < 3; i3++)
280 {
281 i = 12 * i1 + 3 * i2 + i3;
282 for (j2 = 0; j2 < 4; j2++)
283 {
284 if (j2 != i2)
285 {
286 j = 12 * i1 + 3 * j2 + i3;
287 x += freq_[i] * generator_(i, j);
288 }
289 }
290 }
291 }
292 }
293
294 setScale(1 / x);
295
296 if (isScalable())
297 {
298 for (i = 0; i < 36; i++)
299 {
300 eigenValues_[i] /= x;
301 iEigenValues_[i] /= x;
302 }
303 }
304
305 if (!isNonSingular_)
307
308 // and the exchangeability_
309 for (i = 0; i < size_; i++)
310 {
311 for (j = 0; j < size_; j++)
312 {
313 exchangeability_(i, j) = generator_(i, j) / freq_[j];
314 }
315 }
316}
317
319{
320 auto alph = pm.getAlphabet();
322 throw Exception("Need a DNA model");
323
324 // Check that the model is good for YpR, ie transversion rates do
325 // not depend on the origin state
326
327 if ((pm.Qij(0, 1) != pm.Qij(2, 1)) || (pm.Qij(0, 3) != pm.Qij(2, 3))
328 || (pm.Qij(1, 0) != pm.Qij(3, 0)) || (pm.Qij(1, 2) != pm.Qij(3, 2)))
329 throw Exception("Not R/Y Model " + pm.getName());
330}
331
332void YpR::setNamespace(const std::string& prefix)
333{
335
336 // We also need to update the namespace of the nested model:
337 pmodel_->setNamespace(prefix + nestedPrefix_);
338}
339
340// ///////////////////////////////////////////////
341// ///////////////////////////////////////////////
342
343
344/******************************************************************************/
345
347 std::shared_ptr<const RNY> alph,
348 std::unique_ptr<NucleotideSubstitutionModelInterface> pm,
349 double CgT, double TgC,
350 double CaT, double TaC) :
351 AbstractParameterAliasable("YpR_Sym."),
352 YpR(alph, std::move(pm), "YpR_Sym.")
353{
354 addParameter_(new Parameter("YpR_Sym.rCgT", CgT, Parameter::R_PLUS));
355 addParameter_(new Parameter("YpR_Sym.rTgC", TgC, Parameter::R_PLUS));
356 addParameter_(new Parameter("YpR_Sym.rCaT", CaT, Parameter::R_PLUS));
357 addParameter_(new Parameter("YpR_Sym.rTaC", TaC, Parameter::R_PLUS));
358
360}
361
363{
364 double rCgT = getParameterValue("rCgT");
365 double rTgC = getParameterValue("rTgC");
366 double rCaT = getParameterValue("rCaT");
367 double rTaC = getParameterValue("rTaC");
368
369 YpR::updateMatrices_(rCgT, rCgT, rTgC, rTgC, rCaT, rCaT, rTaC, rTaC);
370}
371
374 YpR(ypr)
375{}
376
377/******************************************************************************/
378
379std::string YpR_Sym::getName() const
380{
381 return "YpR_Sym";
382}
383
384
385// ///////////////////////////////////////////////
386// ///////////////////////////////////////////////
387
388/******************************************************************************/
389
391 std::shared_ptr<const RNY> alph,
392 std::unique_ptr<NucleotideSubstitutionModelInterface> pm,
393 double CgT, double cGA,
394 double TgC, double tGA,
395 double CaT, double cAG,
396 double TaC, double tAG) :
397 AbstractParameterAliasable("YpR_Gen."),
398 YpR(alph, std::move(pm), "YpR_Gen.")
399{
400 addParameter_(new Parameter("YpR_Gen.rCgT", CgT, Parameter::R_PLUS));
401 addParameter_(new Parameter("YpR_Gen.rcGA", cGA, Parameter::R_PLUS));
402 addParameter_(new Parameter("YpR_Gen.rTgC", TgC, Parameter::R_PLUS));
403 addParameter_(new Parameter("YpR_Gen.rtGA", tGA, Parameter::R_PLUS));
404 addParameter_(new Parameter("YpR_Gen.rCaT", CaT, Parameter::R_PLUS));
405 addParameter_(new Parameter("YpR_Gen.rcAG", cAG, Parameter::R_PLUS));
406 addParameter_(new Parameter("YpR_Gen.rTaC", TaC, Parameter::R_PLUS));
407 addParameter_(new Parameter("YpR_Gen.rtAG", tAG, Parameter::R_PLUS));
408
410}
411
413{
414 double rCgT = getParameterValue("rCgT");
415 double rcGA = getParameterValue("rcGA");
416 double rTgC = getParameterValue("rTgC");
417 double rtGA = getParameterValue("rtGA");
418 double rCaT = getParameterValue("rCaT");
419 double rcAG = getParameterValue("rcAG");
420 double rTaC = getParameterValue("rTaC");
421 double rtAG = getParameterValue("rtAG");
422
423 YpR::updateMatrices_(rCgT, rcGA, rTgC, rtGA, rCaT, rcAG, rTaC, rtAG);
424}
425
428 YpR(ypr)
429{
431}
432
433/******************************************************************************/
434
435std::string YpR_Gen::getName() const
436{
437 return "YpR_Gen";
438}
void addParameters_(const ParameterList &parameters)
void addParameter_(Parameter *parameter)
void setNamespace(const std::string &prefix)
double getParameterValue(const std::string &name) const override
RowMatrix< double > generator_
The generator matrix of the model.
bool isDiagonalizable_
boolean value for diagonalizability in R of the generator_
bool isNonSingular_
boolean value for non-singularity of rightEigenVectors_
std::vector< RowMatrix< double > > vPowGen_
vector of the powers of generator_ for Taylor development (if rightEigenVectors_ is singular).
Vdouble eigenValues_
The vector of eigen values.
RowMatrix< double > exchangeability_
The exchangeability matrix of the model, defined as . When the model is reversible,...
Vdouble iEigenValues_
The vector of the imaginary part of the eigen values.
RowMatrix< double > leftEigenVectors_
The matrix made of left eigen vectors (by row) if rightEigenVectors_ is non-singular.
virtual bool isScalable() const
returns if model is scalable
RowMatrix< double > rightEigenVectors_
The matrix made of right eigen vectors (by column).
void setScale(double scale)
Multiplies the current generator by the given scale.
RowMatrix< double > tmpMat_
For computational issues.
size_t size_
The number of states.
Vdouble freq_
The vector of equilibrium frequencies.
static bool isNucleicAlphabet(const Alphabet &alphabet)
static void displayMessage(const std::string &text)
virtual std::shared_ptr< const Alphabet > getAlphabet() const =0
virtual std::string getName() const =0
Get the name of the model.
This class implements a state map where all resolved states are modeled.
Definition: StateMap.h:168
const RowMatrix< Real > & getV() const
const std::vector< Real > & getImagEigenValues() const
const std::vector< Real > & getRealEigenValues() const
static Scalar inv(const Matrix< Scalar > &A, Matrix< Scalar > &O)
static void Taylor(const Matrix &A, size_t p, std::vector< RowMatrix< Scalar > > &vO)
static void pow(const Matrix &A, size_t p, Matrix &O)
static void add(MatrixA &A, const MatrixB &B)
static void getId(size_t n, Matrix &O)
static double TINY()
static double SMALL()
static const std::shared_ptr< IntervalConstraint > R_PLUS
Interface for all substitution models.
virtual double Qij(size_t i, size_t j) const =0
A method for computing all necessary matrices.
General YpR model.
Definition: YpR.h:212
void updateMatrices_() override
Diagonalize the matrix, and fill the eigenValues_, iEigenValues_, leftEigenVectors_ and rightEigenVe...
Definition: YpR.cpp:412
std::string getName() const override
Get the name of the model.
Definition: YpR.cpp:435
YpR_Gen(std::shared_ptr< const RNY > alph, std::unique_ptr< NucleotideSubstitutionModelInterface > pm, double CgT=0., double cGA=0., double TgC=0., double tGA=0., double CaT=0., double cAG=0., double TaC=0., double tAG=0.)
Build a new YpR_Gen substitution model.
Definition: YpR.cpp:390
symmetrical YpR model.
Definition: YpR.h:171
std::string getName() const override
Get the name of the model.
Definition: YpR.cpp:379
YpR_Sym(std::shared_ptr< const RNY > alph, std::unique_ptr< NucleotideSubstitutionModelInterface > pm, double CgT=0., double TgC=0., double CaT=0., double TaC=0.)
Build a new YpR_Sym substitution model.
Definition: YpR.cpp:346
void updateMatrices_() override
Diagonalize the matrix, and fill the eigenValues_, iEigenValues_, leftEigenVectors_ and rightEigenVe...
Definition: YpR.cpp:362
YpR model.
Definition: YpR.h:88
std::string nestedPrefix_
Definition: YpR.h:93
void checkModel(const SubstitutionModelInterface &model) const
Definition: YpR.cpp:318
virtual void updateMatrices_() override
Diagonalize the matrix, and fill the eigenValues_, iEigenValues_, leftEigenVectors_ and rightEigenVe...
Definition: YpR.cpp:46
YpR(std::shared_ptr< const RNY >, std::unique_ptr< NucleotideSubstitutionModelInterface > const, const std::string &prefix)
Build a new YpR substitution model, with no dependency parameters.
Definition: YpR.cpp:22
std::unique_ptr< NucleotideSubstitutionModelInterface > pmodel_
Definition: YpR.h:90
virtual void setNamespace(const std::string &) override
Definition: YpR.cpp:332
Defines the basic types of data flow nodes.
ExtendedFloat abs(const ExtendedFloat &ef)