bpp-phyl  2.4.0
bpp::TS98 Class Referenceabstract

Tuffley and Steel's 1998 covarion model. More...

#include <Bpp/Phyl/Model/TS98.h>

+ Inheritance diagram for bpp::TS98:
+ Collaboration diagram for bpp::TS98:

Public Member Functions

 TS98 (ReversibleSubstitutionModel *model, double s1=1., double s2=1., bool normalizeRateChanges=false)
 Build a new TS98 substitution model. More...
 
virtual ~TS98 ()
 
TS98clone () const
 
std::string getName () const
 Get the name of the model. More...
 
double getRate () const
 Get the rate. More...
 
void setRate (double rate)
 Set the rate of the model (must be positive). More...
 
void addRateParameter ()
 
const AlphabetgetAlphabet () const
 
size_t getNumberOfStates () const
 Get the number of states. More...
 
const StateMapgetStateMap () const
 
const std::vector< int > & getAlphabetStates () const
 
std::string getAlphabetStateAsChar (size_t index) const
 
int getAlphabetStateAsInt (size_t index) const
 
std::vector< size_t > getModelStates (int code) const
 Get the state in the model corresponding to a particular state in the alphabet. More...
 
std::vector< size_t > getModelStates (const std::string &code) const
 Get the state in the model corresponding to a particular state in the alphabet. More...
 
const VdoublegetFrequencies () const
 
const Matrix< double > & getExchangeabilityMatrix () const
 
const Matrix< double > & getGenerator () const
 
const Matrix< double > & getPij_t (double t) const
 
const Matrix< double > & getdPij_dt (double t) const
 
const Matrix< double > & getd2Pij_dt2 (double t) const
 
const VdoublegetEigenValues () const
 
const VdoublegetIEigenValues () const
 
bool isDiagonalizable () const
 
bool isNonSingular () const
 
const Matrix< double > & getRowLeftEigenVectors () const
 
const Matrix< double > & getColumnRightEigenVectors () const
 
double freq (size_t i) const
 
double Sij (size_t i, size_t j) const
 
double Qij (size_t i, size_t j) const
 
double Pij_t (size_t i, size_t j, double t) const
 
double dPij_dt (size_t i, size_t j, double t) const
 
double d2Pij_dt2 (size_t i, size_t j, double t) const
 
double getInitValue (size_t i, int state) const
 
void setFreqFromData (const SequenceContainer &data, double pseudoCount=0)
 Set equilibrium frequencies equal to the frequencies estimated from the data. More...
 
virtual void setFreq (std::map< int, double > &frequencies)
 Set equilibrium frequencies. More...
 
const ReversibleSubstitutionModelgetNestedModel () const
 
size_t getRate (size_t i) const
 Get the rate category corresponding to a particular state in the compound model. More...
 
bool isScalable () const
 returns if model is scalable More...
 
void setScalable (bool scalable)
 sets if model is scalable, ie scale can be changed. Default : true, set to false to avoid normalization for example. More...
 
void normalize ()
 Normalize the generator. More...
 
void setDiagonal ()
 set the diagonal of the generator such that sum on each line equals 0. More...
 
double getScale () const
 Get the scalar product of diagonal elements of the generator and the frequencies vector. If the generator is normalized, then scale=1. Otherwise each element must be multiplied by 1/scale. More...
 
void setScale (double scale)
 Multiplies the current generator by the given scale. More...
 
void enableEigenDecomposition (bool yn)
 Set if eigenValues and Vectors must be computed. More...
 
bool enableEigenDecomposition ()
 Tell if eigenValues and Vectors must be computed. More...
 
bool computeFrequencies () const
 
void computeFrequencies (bool yn)
 
virtual void fireParameterChanged (const ParameterList &parameters)
 Tells the model that a parameter value has changed. More...
 
void setNamespace (const std::string &prefix)
 
virtual const FrequenciesSetgetFrequenciesSet () const
 If the model owns a FrequenciesSet, returns a pointer to it, otherwise return 0. More...
 
virtual size_t getNumberOfIndependentParameters () const =0
 
size_t getNumberOfIndependentParameters () const
 
virtual void aliasParameters (const std::string &p1, const std::string &p2)=0
 
virtual void aliasParameters (std::map< std::string, std::string > &unparsedParams, bool verbose)=0
 
void aliasParameters (const std::string &p1, const std::string &p2)
 
void aliasParameters (std::map< std::string, std::string > &unparsedParams, bool verbose)
 
virtual void unaliasParameters (const std::string &p1, const std::string &p2)=0
 
void unaliasParameters (const std::string &p1, const std::string &p2)
 
virtual const ParameterListgetIndependentParameters () const =0
 
const ParameterListgetIndependentParameters () const
 
virtual std::vector< std::string > getAlias (const std::string &name) const =0
 
virtual std::vector< std::string > getAlias (const std::string &name) const
 
virtual std::map< std::string, std::string > getAliases () const =0
 
virtual std::map< std::string, std::string > getAliases () const
 
virtual bool hasParameter (const std::string &name) const =0
 
bool hasParameter (const std::string &name) const
 
virtual const ParameterListgetParameters () const =0
 
const ParameterListgetParameters () const
 
virtual const ParametergetParameter (const std::string &name) const =0
 
const ParametergetParameter (const std::string &name) const
 
virtual double getParameterValue (const std::string &name) const =0
 
double getParameterValue (const std::string &name) const
 
virtual void setAllParametersValues (const ParameterList &parameters)=0
 
void setAllParametersValues (const ParameterList &parameters)
 
virtual void setParameterValue (const std::string &name, double value)=0
 
void setParameterValue (const std::string &name, double value)
 
virtual void setParametersValues (const ParameterList &parameters)=0
 
void setParametersValues (const ParameterList &parameters)
 
virtual bool matchParametersValues (const ParameterList &parameters)=0
 
bool matchParametersValues (const ParameterList &parameters)
 
virtual size_t getNumberOfParameters () const =0
 
size_t getNumberOfParameters () const
 
virtual std::string getNamespace () const =0
 
std::string getNamespace () const
 
virtual std::string getParameterNameWithoutNamespace (const std::string &name) const =0
 
std::string getParameterNameWithoutNamespace (const std::string &name) const
 
bool hasIndependentParameter (const std::string &name) const
 
ParameterList getAliasedParameters (const ParameterList &pl) const
 
ParameterList getFromParameters (const ParameterList &pl) const
 
std::string getFrom (const std::string &name) const
 
const std::shared_ptr< Parameter > & getSharedParameter (const std::string &name) const
 

Protected Member Functions

void updateRatesModel ()
 Update the rates vector, generator and equilibrium frequencies. More...
 
virtual void updateMatrices ()
 A method for computing all necessary matrices. More...
 
VdoublegetFrequencies_ ()
 
const std::shared_ptr< Parameter > & getSharedParameter (size_t i) const
 
std::shared_ptr< Parameter > & getSharedParameter (size_t i)
 
void addParameter_ (Parameter *parameter)
 
void addParameters_ (const ParameterList &parameters)
 
void includeParameters_ (const ParameterList &parameters)
 
void deleteParameter_ (size_t index)
 
void deleteParameter_ (std::string &name)
 
void deleteParameters_ (const std::vector< std::string > &names)
 
void resetParameters_ ()
 
void shareParameter_ (const std::shared_ptr< Parameter > &parameter)
 
void shareParameters_ (const ParameterList &parameters)
 
ParametergetParameter_ (const std::string &name)
 
ParametergetParameter_ (size_t index)
 
const ParametergetParameter_ (size_t index) const
 
ParametergetParameterWithNamespace_ (const std::string &name)
 
const ParametergetParameterWithNamespace_ (const std::string &name) const
 
ParameterListgetParameters_ ()
 

Protected Attributes

ReversibleSubstitutionModelmodel_
 
MarkovModulatedStateMap stateMap_
 
size_t nbStates_
 
size_t nbRates_
 
RowMatrix< double > ratesGenerator_
 
RowMatrix< double > generator_
 The generator matrix $Q$ of the model. More...
 
RowMatrix< double > exchangeability_
 The exchangeability matrix $S$ of the model. More...
 
RowMatrix< double > leftEigenVectors_
 The $U$ matrix made of left eigen vectors (by row). More...
 
RowMatrix< double > rightEigenVectors_
 The $U^-1$ matrix made of right eigen vectors (by column). More...
 
Vdouble eigenValues_
 The vector of real parts of eigen values. More...
 
Vdouble iEigenValues_
 The vector of imaginary parts of the eigen values (zero in case of reversible pmodel). More...
 
bool eigenDecompose_
 Tell if the eigen decomposition should be performed. More...
 
bool compFreq_
 Tell if the equilibrium frequencies should be computed from the generator. More...
 
RowMatrix< double > pijt_
 These ones are for bookkeeping: More...
 
RowMatrix< double > dpijt_
 
RowMatrix< double > d2pijt_
 
Vdouble freq_
 The vector of equilibrium frequencies. More...
 
bool normalizeRateChanges_
 
std::string nestedPrefix_
 
Rate generator.

These variables must be initialized in the constructor of the derived class.

RowMatrix< double > rates_
 
RowMatrix< double > ratesExchangeability_
 
Vdouble ratesFreq_
 

Detailed Description

Tuffley and Steel's 1998 covarion model.

This model is a subclass of the so-called Markov-modulated substitution models, with a rate matrix

\[ G = \begin{pmatrix} -s_1 & s_1\\ s_2 & -s_2 \end{pmatrix} \]

and

\[ D_R = \begin{pmatrix} 0 & 0\\ 0 & \dfrac{s_1+s_2}{s_1} \end{pmatrix}. \]

This model was originally designed for nucleotides sequences, but it can be used with other alphabets.

See also
MarkovModulatedSubstitutionModel

Tuffley C. and Steel M. A., Modelling the covarion hypothesis of nucleotide substitution (1998), Math. Biosci., 147:63-91.

Definition at line 77 of file TS98.h.

Constructor & Destructor Documentation

bpp::TS98::TS98 ( ReversibleSubstitutionModel model,
double  s1 = 1.,
double  s2 = 1.,
bool  normalizeRateChanges = false 
)
inline

Build a new TS98 substitution model.

Parameters
modelThe substitution model to use. May be of any alphabet type.
s1First rate parameter.
s2Second rate parameter.
normalizeRateChangesTell if the rate transition matrix should be normalized.

Definition at line 89 of file TS98.h.

References bpp::AbstractParameterAliasable::addParameter_(), bpp::Parameter::R_PLUS_STAR, bpp::MarkovModulatedSubstitutionModel::updateMatrices(), and updateRatesModel().

Referenced by clone().

virtual bpp::TS98::~TS98 ( )
inlinevirtual

Definition at line 98 of file TS98.h.

Member Function Documentation

void bpp::TS98::addRateParameter ( )
inlinevirtual

Implements bpp::TransitionModel.

Definition at line 109 of file TS98.h.

TS98* bpp::TS98::clone ( ) const
inlinevirtual

Implements bpp::MarkovModulatedSubstitutionModel.

Definition at line 100 of file TS98.h.

References TS98().

bool bpp::MarkovModulatedSubstitutionModel::computeFrequencies ( ) const
inlinevirtualinherited
Returns
Says if equilibrium frequencies should be computed from the generator

Implements bpp::TransitionModel.

Definition at line 306 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::compFreq_.

void bpp::MarkovModulatedSubstitutionModel::computeFrequencies ( bool  yn)
inlinevirtualinherited
Returns
Set if equilibrium frequencies should be computed from the generator

Implements bpp::TransitionModel.

Definition at line 308 of file MarkovModulatedSubstitutionModel.h.

double bpp::MarkovModulatedSubstitutionModel::d2Pij_dt2 ( size_t  i,
size_t  j,
double  t 
) const
inlinevirtualinherited
Returns
The second order derivative of the probability of change from state i to state j with respect to time t, at time t.
See also
getd2Pij_dt2(), getStates()

Implements bpp::TransitionModel.

Definition at line 238 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::getd2Pij_dt2(), and bpp::MarkovModulatedSubstitutionModel::getInitValue().

double bpp::MarkovModulatedSubstitutionModel::dPij_dt ( size_t  i,
size_t  j,
double  t 
) const
inlinevirtualinherited
Returns
The first order derivative of the probability of change from state i to state j with respect to time t, at time t.
See also
getdPij_dt(), getStates()

Implements bpp::TransitionModel.

Definition at line 237 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::getdPij_dt().

void bpp::MarkovModulatedSubstitutionModel::enableEigenDecomposition ( bool  yn)
inlinevirtualinherited

Set if eigenValues and Vectors must be computed.

Implements bpp::SubstitutionModel.

Definition at line 302 of file MarkovModulatedSubstitutionModel.h.

bool bpp::MarkovModulatedSubstitutionModel::enableEigenDecomposition ( )
inlinevirtualinherited

Tell if eigenValues and Vectors must be computed.

Implements bpp::SubstitutionModel.

Definition at line 304 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::eigenDecompose_.

virtual void bpp::MarkovModulatedSubstitutionModel::fireParameterChanged ( const ParameterList parameters)
inlinevirtualinherited
double bpp::MarkovModulatedSubstitutionModel::freq ( size_t  i) const
inlinevirtualinherited
Returns
Equilibrium frequency associated to character i.
See also
getFrequencies(), getStates()

Implements bpp::TransitionModel.

Definition at line 232 of file MarkovModulatedSubstitutionModel.h.

const Alphabet* bpp::MarkovModulatedSubstitutionModel::getAlphabet ( ) const
inlinevirtualinherited
Returns
Get the alphabet associated to this model.

Implements bpp::TransitionModel.

Definition at line 197 of file MarkovModulatedSubstitutionModel.h.

References bpp::TransitionModel::getAlphabet().

Referenced by bpp::MarkovModulatedSubstitutionModel::getInitValue().

std::string bpp::MarkovModulatedSubstitutionModel::getAlphabetStateAsChar ( size_t  index) const
inlinevirtualinherited
Parameters
indexThe model state.
Returns
The corresponding alphabet state as character code.

Implements bpp::TransitionModel.

Definition at line 205 of file MarkovModulatedSubstitutionModel.h.

References bpp::AbstractStateMap::getAlphabetStateAsChar().

int bpp::MarkovModulatedSubstitutionModel::getAlphabetStateAsInt ( size_t  index) const
inlinevirtualinherited
Parameters
indexThe model state.
Returns
The corresponding alphabet state as character code.

Implements bpp::TransitionModel.

Definition at line 207 of file MarkovModulatedSubstitutionModel.h.

References bpp::AbstractStateMap::getAlphabetStateAsInt().

Referenced by bpp::MarkovModulatedSubstitutionModel::getInitValue().

const std::vector<int>& bpp::MarkovModulatedSubstitutionModel::getAlphabetStates ( ) const
inlinevirtualinherited
Returns
The alphabet states of each state of the model, as a vector of int codes.
See also
Alphabet

Implements bpp::TransitionModel.

Definition at line 203 of file MarkovModulatedSubstitutionModel.h.

References bpp::AbstractStateMap::getAlphabetStates().

const Matrix<double>& bpp::MarkovModulatedSubstitutionModel::getColumnRightEigenVectors ( ) const
inlinevirtualinherited
Returns
A matrix with right eigen vectors. Each column in the matrix stands for an eigen vector.

Implements bpp::SubstitutionModel.

Definition at line 230 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::rightEigenVectors_.

const Matrix< double > & MarkovModulatedSubstitutionModel::getd2Pij_dt2 ( double  t) const
virtualinherited
const Matrix< double > & MarkovModulatedSubstitutionModel::getdPij_dt ( double  t) const
virtualinherited
const Vdouble& bpp::MarkovModulatedSubstitutionModel::getEigenValues ( ) const
inlinevirtualinherited
Returns
A vector with all real parts of the eigen values of the generator of this model;

Implements bpp::SubstitutionModel.

Definition at line 223 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::eigenValues_.

const Matrix<double>& bpp::MarkovModulatedSubstitutionModel::getExchangeabilityMatrix ( ) const
inlinevirtualinherited
Returns
The matrix of exchangeability terms. It is recommended that exchangeability matrix be normalized so that the normalized generator be obtained directly by the dot product $S . \pi$.

Implements bpp::SubstitutionModel.

Definition at line 215 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::exchangeability_.

const Vdouble& bpp::MarkovModulatedSubstitutionModel::getFrequencies ( ) const
inlinevirtualinherited
Returns
A vector of all equilibrium frequencies.
See also
freq()

Implements bpp::TransitionModel.

Definition at line 213 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::freq_.

Vdouble& bpp::MarkovModulatedSubstitutionModel::getFrequencies_ ( )
inlineprotectedvirtualinherited
const Matrix<double>& bpp::MarkovModulatedSubstitutionModel::getGenerator ( ) const
inlinevirtualinherited
Returns
The normalized Markov generator matrix, i.e. all normalized rates of changes from state i to state j. The generator is normalized so that (i) $ \forall i; \sum_j Q_{i,j} = 0 $, meaning that $ $ \forall i; Q_{i,i} = -\sum_{j \neq i}Q_{i,j}$, and (ii) $ \sum_i Q_{i,i} \times \pi_i = -1$. This means that, under normalization, the mean rate of replacement at equilibrium is 1 and that time $t$ are measured in units of expected number of changes per site. Additionnaly, the rate_ attibute provides the possibility to increase or decrease this mean rate.

See Kosiol and Goldman (2005), Molecular Biology And Evolution 22(2) 193-9.

See also
Qij()

Implements bpp::SubstitutionModel.

Definition at line 217 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::generator_, bpp::MarkovModulatedSubstitutionModel::getd2Pij_dt2(), bpp::MarkovModulatedSubstitutionModel::getdPij_dt(), and bpp::MarkovModulatedSubstitutionModel::getPij_t().

const Vdouble& bpp::MarkovModulatedSubstitutionModel::getIEigenValues ( ) const
inlinevirtualinherited
Returns
A vector with all imaginary parts of the eigen values of the generator of this model;

Implements bpp::SubstitutionModel.

Definition at line 224 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::iEigenValues_.

double MarkovModulatedSubstitutionModel::getInitValue ( size_t  i,
int  state 
) const
virtualinherited

This method is used to initialize likelihoods in reccursions. It typically sends 1 if i = state, 0 otherwise, where i is one of the possible states of the alphabet allowed in the model and state is the observed state in the considered sequence/site.

Parameters
ithe index of the state in the model.
stateAn observed state in the sequence/site.
Returns
1 or 0 depending if the two states are compatible.
Exceptions
IndexOutOfBoundsExceptionif array position is out of range.
BadIntExceptionif states are not allowed in the associated alphabet.
See also
getStates();

Implements bpp::TransitionModel.

Definition at line 220 of file MarkovModulatedSubstitutionModel.cpp.

References bpp::Alphabet::getAlias(), bpp::MarkovModulatedSubstitutionModel::getAlphabet(), bpp::TransitionModel::getAlphabet(), bpp::MarkovModulatedSubstitutionModel::getAlphabetStateAsInt(), bpp::Alphabet::intToChar(), bpp::MarkovModulatedSubstitutionModel::model_, bpp::MarkovModulatedSubstitutionModel::nbRates_, and bpp::MarkovModulatedSubstitutionModel::nbStates_.

Referenced by bpp::MarkovModulatedSubstitutionModel::d2Pij_dt2().

std::vector<size_t> bpp::MarkovModulatedSubstitutionModel::getModelStates ( int  code) const
inlinevirtualinherited

Get the state in the model corresponding to a particular state in the alphabet.

Parameters
codeThe alphabet state to check.
Returns
A vector of indices of model states.

Implements bpp::TransitionModel.

Definition at line 209 of file MarkovModulatedSubstitutionModel.h.

References bpp::AbstractStateMap::getModelStates().

std::vector<size_t> bpp::MarkovModulatedSubstitutionModel::getModelStates ( const std::string &  code) const
inlinevirtualinherited

Get the state in the model corresponding to a particular state in the alphabet.

Parameters
codeThe alphabet state to check.
Returns
A vector of indices of model states.

Implements bpp::TransitionModel.

Definition at line 211 of file MarkovModulatedSubstitutionModel.h.

References bpp::AbstractStateMap::getModelStates().

std::string bpp::TS98::getName ( ) const
inlinevirtual

Get the name of the model.

Returns
The name of this model.

Implements bpp::TransitionModel.

Definition at line 103 of file TS98.h.

const ReversibleSubstitutionModel* bpp::MarkovModulatedSubstitutionModel::getNestedModel ( ) const
inlineinherited
size_t bpp::MarkovModulatedSubstitutionModel::getNumberOfStates ( ) const
inlinevirtualinherited

Get the number of states.

For most models, this equals the size of the alphabet.

See also
getAlphabetChars for the list of supported states.
Returns
The number of different states in the model.

Implements bpp::TransitionModel.

Definition at line 199 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::nbRates_.

Referenced by bpp::MarkovModulatedSubstitutionModel::setDiagonal().

double bpp::TS98::getRate ( ) const
inlinevirtual

Get the rate.

Reimplemented from bpp::MarkovModulatedSubstitutionModel.

Definition at line 105 of file TS98.h.

size_t bpp::MarkovModulatedSubstitutionModel::getRate ( size_t  i) const
inlineinherited

Get the rate category corresponding to a particular state in the compound model.

Parameters
iThe state.
Returns
The corresponding rate category.
See also
getState;

Definition at line 263 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::nbStates_.

const Matrix<double>& bpp::MarkovModulatedSubstitutionModel::getRowLeftEigenVectors ( ) const
inlinevirtualinherited
Returns
A matrix with left eigen vectors. Each row in the matrix stands for an eigen vector.

Implements bpp::SubstitutionModel.

Definition at line 229 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::leftEigenVectors_.

double bpp::MarkovModulatedSubstitutionModel::getScale ( ) const
inlinevirtualinherited

Get the scalar product of diagonal elements of the generator and the frequencies vector. If the generator is normalized, then scale=1. Otherwise each element must be multiplied by 1/scale.

Returns
Minus the scalar product of diagonal elements and the frequencies vector.

Implements bpp::SubstitutionModel.

Definition at line 290 of file MarkovModulatedSubstitutionModel.h.

References bpp::MatrixTools::diag(), and bpp::MarkovModulatedSubstitutionModel::freq_.

Referenced by bpp::MarkovModulatedSubstitutionModel::updateMatrices().

const StateMap& bpp::MarkovModulatedSubstitutionModel::getStateMap ( ) const
inlinevirtualinherited
Returns
The mapping of model states with alphabet states.

Implements bpp::TransitionModel.

Definition at line 201 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::stateMap_.

bool bpp::MarkovModulatedSubstitutionModel::isDiagonalizable ( ) const
inlinevirtualinherited
Returns
True if the model is diagonalizable in R.

Implements bpp::SubstitutionModel.

Definition at line 226 of file MarkovModulatedSubstitutionModel.h.

bool bpp::MarkovModulatedSubstitutionModel::isNonSingular ( ) const
inlinevirtualinherited
Returns
True is the model is non-singular.

Implements bpp::SubstitutionModel.

Definition at line 227 of file MarkovModulatedSubstitutionModel.h.

bool bpp::MarkovModulatedSubstitutionModel::isScalable ( ) const
inlinevirtualinherited

returns if model is scalable

Implements bpp::SubstitutionModel.

Definition at line 272 of file MarkovModulatedSubstitutionModel.h.

References bpp::SubstitutionModel::isScalable().

Referenced by bpp::MarkovModulatedSubstitutionModel::updateMatrices().

void bpp::MarkovModulatedSubstitutionModel::normalize ( )
inlinevirtualinherited
double bpp::MarkovModulatedSubstitutionModel::Pij_t ( size_t  i,
size_t  j,
double  t 
) const
inlinevirtualinherited
Returns
The probability of change from state i to state j during time t.
See also
getPij_t(), getStates()

Implements bpp::TransitionModel.

Definition at line 236 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::getPij_t().

double bpp::MarkovModulatedSubstitutionModel::Qij ( size_t  i,
size_t  j 
) const
inlinevirtualinherited
Returns
The rate in the generator of change from state i to state j.
See also
getStates();

Implements bpp::SubstitutionModel.

Definition at line 234 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::generator_.

void MarkovModulatedSubstitutionModel::setDiagonal ( )
virtualinherited
virtual void bpp::MarkovModulatedSubstitutionModel::setFreq ( std::map< int, double > &  frequencies)
inlinevirtualinherited

Set equilibrium frequencies.

Parameters
frequenciesThe map of the frequencies to use.

Implements bpp::TransitionModel.

Definition at line 248 of file MarkovModulatedSubstitutionModel.h.

References bpp::TransitionModel::setFreq(), and bpp::MarkovModulatedSubstitutionModel::updateMatrices().

void bpp::MarkovModulatedSubstitutionModel::setFreqFromData ( const SequenceContainer data,
double  pseudoCount = 0 
)
inlinevirtualinherited

Set equilibrium frequencies equal to the frequencies estimated from the data.

Parameters
dataThe sequences to use.
pseudoCountA quantity $\psi$ to add to adjust the observed values in order to prevent issues due to missing states on small data set. The corrected frequencies shall be computed as

\[ \pi_i = \frac{n_i+\psi}{\sum_j (f_j+\psi)} \]

Implements bpp::TransitionModel.

Definition at line 242 of file MarkovModulatedSubstitutionModel.h.

References bpp::TransitionModel::setFreqFromData(), and bpp::MarkovModulatedSubstitutionModel::updateMatrices().

void bpp::TS98::setRate ( double  rate)
inlinevirtual

Set the rate of the model (must be positive).

Parameters
ratemust be positive.

Reimplemented from bpp::MarkovModulatedSubstitutionModel.

Definition at line 107 of file TS98.h.

void bpp::MarkovModulatedSubstitutionModel::setScalable ( bool  scalable)
inlinevirtualinherited

sets if model is scalable, ie scale can be changed. Default : true, set to false to avoid normalization for example.

Implements bpp::SubstitutionModel.

Definition at line 277 of file MarkovModulatedSubstitutionModel.h.

References bpp::SubstitutionModel::setScalable().

void bpp::MarkovModulatedSubstitutionModel::setScale ( double  scale)
inlinevirtualinherited

Multiplies the current generator by the given scale.

Parameters
scalethe scale by which the generator is multiplied.

Implements bpp::SubstitutionModel.

Definition at line 297 of file MarkovModulatedSubstitutionModel.h.

References bpp::SubstitutionModel::setScale(), and bpp::MarkovModulatedSubstitutionModel::updateMatrices().

Referenced by bpp::MarkovModulatedSubstitutionModel::updateMatrices().

double bpp::MarkovModulatedSubstitutionModel::Sij ( size_t  i,
size_t  j 
) const
inlinevirtualinherited
Returns
The exchangeability between state i and state j.

By definition Sij(i,j) = Sij(j,i).

Implements bpp::SubstitutionModel.

Definition at line 233 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::exchangeability_.

void MarkovModulatedSubstitutionModel::updateMatrices ( )
protectedvirtualinherited

A method for computing all necessary matrices.

Implements bpp::SubstitutionModel.

Definition at line 109 of file MarkovModulatedSubstitutionModel.cpp.

References bpp::MatrixTools::add(), bpp::MarkovModulatedSubstitutionModel::d2pijt_, bpp::MatrixTools::diag(), bpp::MarkovModulatedSubstitutionModel::dpijt_, bpp::MarkovModulatedSubstitutionModel::eigenValues_, bpp::MarkovModulatedSubstitutionModel::exchangeability_, bpp::MarkovModulatedSubstitutionModel::freq_, bpp::MarkovModulatedSubstitutionModel::generator_, bpp::SubstitutionModel::getColumnRightEigenVectors(), bpp::SubstitutionModel::getEigenValues(), bpp::SubstitutionModel::getExchangeabilityMatrix(), bpp::TransitionModel::getFrequencies(), bpp::SubstitutionModel::getGenerator(), bpp::RowMatrix< Scalar >::getNumberOfColumns(), bpp::TransitionModel::getNumberOfStates(), bpp::EigenValue< Real >::getRealEigenValues(), bpp::MarkovModulatedSubstitutionModel::getScale(), bpp::EigenValue< Real >::getV(), bpp::MarkovModulatedSubstitutionModel::iEigenValues_, bpp::MatrixTools::inv(), bpp::MarkovModulatedSubstitutionModel::isScalable(), bpp::VectorTools::kroneckerMult(), bpp::MatrixTools::kroneckerMult(), bpp::MarkovModulatedSubstitutionModel::leftEigenVectors_, bpp::MarkovModulatedSubstitutionModel::model_, bpp::MatrixTools::mult(), bpp::MarkovModulatedSubstitutionModel::nbRates_, bpp::MarkovModulatedSubstitutionModel::nbStates_, bpp::MarkovModulatedSubstitutionModel::normalizeRateChanges_, bpp::MarkovModulatedSubstitutionModel::pijt_, bpp::MarkovModulatedSubstitutionModel::rates_, bpp::MarkovModulatedSubstitutionModel::ratesExchangeability_, bpp::MarkovModulatedSubstitutionModel::ratesFreq_, bpp::MarkovModulatedSubstitutionModel::ratesGenerator_, bpp::RowMatrix< Scalar >::resize(), bpp::MarkovModulatedSubstitutionModel::rightEigenVectors_, bpp::MatrixTools::scale(), and bpp::MarkovModulatedSubstitutionModel::setScale().

Referenced by bpp::MarkovModulatedSubstitutionModel::fireParameterChanged(), bpp::G2001::G2001(), bpp::MarkovModulatedSubstitutionModel::normalize(), bpp::MarkovModulatedSubstitutionModel::setFreq(), bpp::MarkovModulatedSubstitutionModel::setFreqFromData(), bpp::MarkovModulatedSubstitutionModel::setScale(), and TS98().

void bpp::TS98::updateRatesModel ( )
inlineprotectedvirtual

Update the rates vector, generator and equilibrium frequencies.

This method must be implemented by the derived class. It is called by the fireParameterChanged() method.

Implements bpp::MarkovModulatedSubstitutionModel.

Definition at line 112 of file TS98.h.

References bpp::AbstractParameterAliasable::getParameterValue(), bpp::MarkovModulatedSubstitutionModel::rates_, bpp::MarkovModulatedSubstitutionModel::ratesExchangeability_, and bpp::MarkovModulatedSubstitutionModel::ratesFreq_.

Referenced by TS98().

Member Data Documentation

bool bpp::MarkovModulatedSubstitutionModel::compFreq_
protectedinherited

Tell if the equilibrium frequencies should be computed from the generator.

Definition at line 145 of file MarkovModulatedSubstitutionModel.h.

Referenced by bpp::MarkovModulatedSubstitutionModel::computeFrequencies(), and bpp::MarkovModulatedSubstitutionModel::operator=().

RowMatrix<double> bpp::MarkovModulatedSubstitutionModel::d2pijt_
mutableprotectedinherited
RowMatrix<double> bpp::MarkovModulatedSubstitutionModel::dpijt_
mutableprotectedinherited
bool bpp::MarkovModulatedSubstitutionModel::eigenDecompose_
protectedinherited

Tell if the eigen decomposition should be performed.

Definition at line 138 of file MarkovModulatedSubstitutionModel.h.

Referenced by bpp::MarkovModulatedSubstitutionModel::enableEigenDecomposition(), and bpp::MarkovModulatedSubstitutionModel::operator=().

RowMatrix<double> bpp::MarkovModulatedSubstitutionModel::exchangeability_
protectedinherited
Vdouble bpp::MarkovModulatedSubstitutionModel::iEigenValues_
protectedinherited

The vector of imaginary parts of the eigen values (zero in case of reversible pmodel).

Definition at line 133 of file MarkovModulatedSubstitutionModel.h.

Referenced by bpp::MarkovModulatedSubstitutionModel::getIEigenValues(), bpp::MarkovModulatedSubstitutionModel::operator=(), and bpp::MarkovModulatedSubstitutionModel::updateMatrices().

std::string bpp::MarkovModulatedSubstitutionModel::nestedPrefix_
protectedinherited
bool bpp::MarkovModulatedSubstitutionModel::normalizeRateChanges_
protectedinherited
RowMatrix<double> bpp::MarkovModulatedSubstitutionModel::pijt_
mutableprotectedinherited
RowMatrix<double> bpp::MarkovModulatedSubstitutionModel::ratesExchangeability_
protectedinherited
Vdouble bpp::MarkovModulatedSubstitutionModel::ratesFreq_
protectedinherited
RowMatrix<double> bpp::MarkovModulatedSubstitutionModel::ratesGenerator_
protectedinherited
MarkovModulatedStateMap bpp::MarkovModulatedSubstitutionModel::stateMap_
protectedinherited

The documentation for this class was generated from the following file: