42 #include "Bpp/Numeric/Matrix/MatrixTools.h" 54 nbStates_(model->getNumberOfStates()),
55 jMat_(nbStates_, nbStates_),
56 v_(nbStates_, nbStates_),
57 vInv_(nbStates_, nbStates_),
58 lambda_(nbStates_, nbStates_),
59 bMatrices_(reg->getNumberOfSubstitutionTypes()),
60 insideProducts_(reg->getNumberOfSubstitutionTypes()),
61 counts_(reg->getNumberOfSubstitutionTypes()),
66 throw Exception(
"DecompositionSubstitutionCount (constructor): alphabets do not match between register and model.");
84 if (i > 0 && k != j) {
102 for (
size_t i = 0; i <
register_->getNumberOfSubstitutionTypes(); ++i) {
112 size_t nbTypes =
register_->getNumberOfSubstitutionTypes();
129 for (
size_t i = 0; i <
register_->getNumberOfSubstitutionTypes(); ++i) {
138 vector<double> expLam = VectorTools::exp(lambda * t);
139 for (
unsigned int i = 0; i <
nbStates_; ++i) {
140 for (
unsigned int j = 0; j <
nbStates_; ++j) {
141 double dd = lambda[i] - lambda[j];
143 result(i, j) = t * expLam[i];
145 result(i, j) = (expLam[i] - expLam[j]) / dd;
156 for (
size_t i = 0; i <
register_->getNumberOfSubstitutionTypes(); ++i) {
159 MatrixTools::mult(
v_, tmp1, tmp2);
165 for (
size_t i = 0; i <
register_->getNumberOfSubstitutionTypes(); i++) {
173 counts_[i](j, k) *=
weights_->getIndex(supportedStates[j], supportedStates[k]);
185 throw Exception(
"DecompositionSubstitutionCount::getAllNumbersOfSubstitutions. Negative branch length: " + TextTools::toString(length) +
".");
207 return new RowMatrix<double>(
counts_[type - 1]);
215 throw Exception(
"DecompositionSubstitutionCount::getNumbersOfSubstitutions. Negative branch length: " + TextTools::toString(length) +
".");
221 return counts_[type - 1](initialState, finalState);
229 throw Exception(
"DecompositionSubstitutionCount::getNumbersOfSubstitutions. Negative branch length: " + TextTools::toString(length) +
".");
237 v[t] =
counts_[t](initialState, finalState);
248 throw Exception(
"DecompositionSubstitutionCount::setSubstitutionModel. Only works with reversible models for now.");
252 throw Exception(
"DecompositionSubstitutionCount::setSubstitutionModel: alphabets do not match between register and model.");
274 throw Exception(
"DecompositionSubstitutionCount::substitutionRegisterHasChanged: alphabets do not match between register and model.");
290 if (
weights_->getAlphabet()->getAlphabetType() !=
register_->getAlphabet()->getAlphabetType())
291 throw Exception(
"DecompositionSubstitutionCount::weightsHaveChanged. Incorrect alphabet type.");
std::vector< double > getNumberOfSubstitutionsForEachType(size_t initialState, size_t finalState, double length) const
Get the numbers of susbstitutions on a branch for all types, for an initial and final states...
Interface for all substitution models.
std::auto_ptr< SubstitutionRegister > register_
const ReversibleSubstitutionModel * model_
RowMatrix< double > vInv_
virtual const std::vector< int > & getAlphabetStates() const =0
Basic implementation of the the SubstitutionCount interface.
double getNumberOfSubstitutions(size_t initialState, size_t finalState, double length, size_t type=1) const
Get the number of susbstitutions on a branch, given the initial and final states, and the branch leng...
RowMatrix< double > jMat_
The SubstitutionRegister interface.
virtual const Alphabet * getAlphabet() const =0
void computeCounts_(double length) const
const AlphabetIndex2 * weights_
void substitutionRegisterHasChanged()
virtual const Vdouble & getEigenValues() const =0
virtual const Matrix< double > & getPij_t(double t) const =0
std::vector< RowMatrix< double > > counts_
Interface for reversible substitution models.
virtual size_t getNumberOfSubstitutionTypes() const
Short cut function, equivalent to getSubstitutionRegister().getNumberOfSubstitutionTypes().
void weightsHaveChanged()
virtual double Qij(size_t i, size_t j) const =0
Matrix< double > * getAllNumbersOfSubstitutions(double length, size_t type=1) const
Get the numbers of susbstitutions on a branch, for each initial and final states, and given the branc...
void setSubstitutionModel(const SubstitutionModel *model)
Set the substitution model.
std::vector< double > lambda_
DecompositionSubstitutionCount(const ReversibleSubstitutionModel *model, SubstitutionRegister *reg, const AlphabetIndex2 *weights=0)
virtual const Matrix< double > & getRowLeftEigenVectors() const =0
virtual const Matrix< double > & getColumnRightEigenVectors() const =0
std::vector< RowMatrix< double > > insideProducts_
Partial implementation of the WeightedSubstitutionCount interface.
void jFunction_(const std::vector< double > &lambda, double t, RowMatrix< double > &result) const
virtual const Alphabet * getAlphabet() const =0
std::vector< RowMatrix< double > > bMatrices_