bpp-phyl3  3.0.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 (std::unique_ptr< ReversibleSubstitutionModelInterface > model, double s1=1., double s2=1., bool normalizeRateChanges=false)
 Build a new TS98 substitution model. More...
 
virtual ~TS98 ()
 
TS98clone () const override
 
std::string getName () const override
 Get the name of the model. More...
 
double getRate () const override
 Get the rate. More...
 
void setRate (double rate) override
 Set the rate of the model (must be positive). More...
 
void addRateParameter () override
 
const Alphabetalphabet () const override
 
std::shared_ptr< const AlphabetgetAlphabet () const override
 
size_t getNumberOfStates () const override
 Get the number of states. More...
 
const StateMapInterfacestateMap () const override
 
std::shared_ptr< const StateMapInterfacegetStateMap () const override
 
const std::vector< int > & getAlphabetStates () const override
 
std::string getAlphabetStateAsChar (size_t index) const override
 
int getAlphabetStateAsInt (size_t index) const override
 
std::vector< size_t > getModelStates (int code) const override
 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 override
 Get the state in the model corresponding to a particular state in the alphabet. More...
 
const VdoublegetFrequencies () const override
 
const Matrix< double > & exchangeabilityMatrix () const override
 
const Matrix< double > & generator () const override
 
const Matrix< double > & getPij_t (double t) const override
 
const Matrix< double > & getdPij_dt (double t) const override
 
const Matrix< double > & getd2Pij_dt2 (double t) const override
 
const VdoublegetEigenValues () const override
 
const VdoublegetIEigenValues () const override
 
bool isDiagonalizable () const override
 
bool isNonSingular () const override
 
const Matrix< double > & getRowLeftEigenVectors () const override
 
const Matrix< double > & getColumnRightEigenVectors () const override
 
double freq (size_t i) const override
 
double Sij (size_t i, size_t j) const override
 
double Qij (size_t i, size_t j) const override
 A method for computing all necessary matrices. More...
 
double Pij_t (size_t i, size_t j, double t) const override
 
double dPij_dt (size_t i, size_t j, double t) const override
 
double d2Pij_dt2 (size_t i, size_t j, double t) const override
 
double getInitValue (size_t i, int state) const override
 
void setFreqFromData (const SequenceDataInterface &data, double pseudoCount=0) override
 Set equilibrium frequencies equal to the frequencies estimated from the data. More...
 
void setFreq (std::map< int, double > &frequencies) override
 Set equilibrium frequencies. More...
 
const FrequencySetInterfacefrequencySet () const override
 
const ReversibleSubstitutionModelInterfacenestedModel () 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 override
 returns if model is scalable More...
 
void setScalable (bool scalable) override
 sets if model is scalable, ie scale can be changed. Default : true, set to false to avoid normalization for example. More...
 
void normalize () override
 Normalize the generator. More...
 
void setDiagonal () override
 set the diagonal of the generator such that sum on each line equals 0. More...
 
double getScale () const override
 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) override
 Multiplies the current generator by the given scale. More...
 
void enableEigenDecomposition (bool yn) override
 Set if eigenValues and Vectors must be computed. More...
 
bool enableEigenDecomposition () override
 Tell if eigenValues and Vectors must be computed. More...
 
bool computeFrequencies () const override
 
void computeFrequencies (bool yn) override
 
virtual void fireParameterChanged (const ParameterList &parameters) override
 Tells the model that a parameter value has changed. More...
 
void setNamespace (const std::string &prefix) override
 
virtual const Eigen::VectorXd & Lik_t (const Eigen::VectorXd &values, double t) const =0
 
const Eigen::VectorXd & Lik_t (const Eigen::VectorXd &values, double t) const override
 
virtual const Eigen::VectorXd & dLik_dt (const Eigen::VectorXd &values, double t) const =0
 
const Eigen::VectorXd & dLik_dt (const Eigen::VectorXd &values, double t) const override
 
virtual const Eigen::VectorXd & d2Lik_dt2 (const Eigen::VectorXd &values, double t) const =0
 
const Eigen::VectorXd & d2Lik_dt2 (const Eigen::VectorXd &values, double t) const override
 
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 override
 
virtual const ParameterListgetParameters () const=0
 
const ParameterListgetParameters () const override
 
virtual const Parameterparameter (const std::string &name) const=0
 
const Parameterparameter (const std::string &name) const override
 
virtual double getParameterValue (const std::string &name) const=0
 
double getParameterValue (const std::string &name) const override
 
virtual void setAllParametersValues (const ParameterList &parameters)=0
 
void setAllParametersValues (const ParameterList &parameters) override
 
virtual void setParameterValue (const std::string &name, double value)=0
 
void setParameterValue (const std::string &name, double value) override
 
virtual void setParametersValues (const ParameterList &parameters)=0
 
void setParametersValues (const ParameterList &parameters) override
 
virtual bool matchParametersValues (const ParameterList &parameters)=0
 
bool matchParametersValues (const ParameterList &parameters) override
 
virtual void removeConstraint (const std::string &name)=0
 
void removeConstraint (const std::string &name) override
 
virtual void setConstraint (const std::string &name, std::shared_ptr< ConstraintInterface > constraint)=0
 
void setConstraint (const std::string &name, std::shared_ptr< ConstraintInterface > constraint) override
 
virtual size_t getNumberOfParameters () const=0
 
size_t getNumberOfParameters () const override
 
virtual std::string getNamespace () const=0
 
std::string getNamespace () const override
 
virtual std::string getParameterNameWithoutNamespace (const std::string &name) const=0
 
std::string getParameterNameWithoutNamespace (const std::string &name) const override
 
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 > & getParameter (const std::string &name) const
 

Protected Member Functions

void updateRatesModel_ () override
 Update the rates vector, generator and equilibrium frequencies. More...
 
virtual void updateMatrices_ ()
 
VdoublegetFrequencies_ () override
 
virtual ParameterListgetParameters_ ()=0
 
ParameterListgetParameters_ () override
 
const std::shared_ptr< Parameter > & getParameter (size_t i) const
 
std::shared_ptr< Parameter > & getParameter (size_t i)
 
void addParameter_ (Parameter *parameter)
 
void addParameters_ (const ParameterList &parameters)
 
void shareParameter_ (const std::shared_ptr< Parameter > &parameter)
 
void shareParameters_ (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_ ()
 
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
 

Protected Attributes

std::unique_ptr< ReversibleSubstitutionModelInterfacemodel_
 
std::shared_ptr< const MarkovModulatedStateMapstateMap_
 
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_
 

Private Attributes

ParameterList independentParameters_
 
std::map< std::string, std::shared_ptr< AliasParameterListener > > aliasListenersRegister_
 
ParameterList parameters_
 
std::string prefix_
 
Eigen::VectorXd lik_
 

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 42 of file TS98.h.

Constructor & Destructor Documentation

◆ TS98()

bpp::TS98::TS98 ( std::unique_ptr< ReversibleSubstitutionModelInterface 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 54 of file TS98.h.

References bpp::AbstractParameterAliasable::addParameter_(), bpp::Parameter::R_PLUS_STAR, bpp::MarkovModulatedSubstitutionModel::updateMatrices_(), and updateRatesModel_().

Referenced by clone().

◆ ~TS98()

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

Definition at line 67 of file TS98.h.

Member Function Documentation

◆ addRateParameter()

void bpp::TS98::addRateParameter ( )
inlineoverridevirtual

Implements bpp::BranchModelInterface.

Definition at line 78 of file TS98.h.

◆ alphabet()

const Alphabet& bpp::MarkovModulatedSubstitutionModel::alphabet ( ) const
inlineoverridevirtualinherited
Returns
A reference to the alphabet associated to this model.

Implements bpp::BranchModelInterface.

Definition at line 177 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::model_.

◆ clone()

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

Implements bpp::MarkovModulatedSubstitutionModel.

Definition at line 69 of file TS98.h.

References TS98().

◆ computeFrequencies() [1/2]

bool bpp::MarkovModulatedSubstitutionModel::computeFrequencies ( ) const
inlineoverridevirtualinherited
Returns
Says if equilibrium frequencies should be computed

Implements bpp::TransitionModelInterface.

Definition at line 299 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::compFreq_.

◆ computeFrequencies() [2/2]

void bpp::MarkovModulatedSubstitutionModel::computeFrequencies ( bool  yn)
inlineoverridevirtualinherited
Returns
Set if equilibrium frequencies should be computed

Implements bpp::TransitionModelInterface.

Definition at line 301 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::compFreq_.

◆ d2Lik_dt2() [1/2]

virtual const Eigen::VectorXd& bpp::TransitionModelInterface::d2Lik_dt2 ( const Eigen::VectorXd &  values,
double  t 
) const
pure virtualinherited

◆ d2Lik_dt2() [2/2]

const Eigen::VectorXd& bpp::AbstractLkTransitionModel::d2Lik_dt2 ( const Eigen::VectorXd &  values,
double  t 
) const
inlineoverridevirtualinherited

◆ d2Pij_dt2()

double bpp::MarkovModulatedSubstitutionModel::d2Pij_dt2 ( size_t  i,
size_t  j,
double  t 
) const
inlineoverridevirtualinherited
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::TransitionModelInterface.

Definition at line 222 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::getd2Pij_dt2().

◆ dLik_dt() [1/2]

virtual const Eigen::VectorXd& bpp::TransitionModelInterface::dLik_dt ( const Eigen::VectorXd &  values,
double  t 
) const
pure virtualinherited

◆ dLik_dt() [2/2]

const Eigen::VectorXd& bpp::AbstractLkTransitionModel::dLik_dt ( const Eigen::VectorXd &  values,
double  t 
) const
inlineoverridevirtualinherited

◆ dPij_dt()

double bpp::MarkovModulatedSubstitutionModel::dPij_dt ( size_t  i,
size_t  j,
double  t 
) const
inlineoverridevirtualinherited
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::TransitionModelInterface.

Definition at line 221 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::getdPij_dt().

◆ enableEigenDecomposition() [1/2]

bool bpp::MarkovModulatedSubstitutionModel::enableEigenDecomposition ( )
inlineoverridevirtualinherited

Tell if eigenValues and Vectors must be computed.

Implements bpp::SubstitutionModelInterface.

Definition at line 297 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::eigenDecompose_.

◆ enableEigenDecomposition() [2/2]

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

Set if eigenValues and Vectors must be computed.

Implements bpp::SubstitutionModelInterface.

Definition at line 295 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::eigenDecompose_.

◆ exchangeabilityMatrix()

const Matrix<double>& bpp::MarkovModulatedSubstitutionModel::exchangeabilityMatrix ( ) const
inlineoverridevirtualinherited
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::SubstitutionModelInterface.

Definition at line 199 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::exchangeability_.

◆ fireParameterChanged()

virtual void bpp::MarkovModulatedSubstitutionModel::fireParameterChanged ( const ParameterList parameters)
inlineoverridevirtualinherited

Tells the model that a parameter value has changed.

This updates the matrices consequently.

Reimplemented from bpp::AbstractParameterAliasable.

Reimplemented in bpp::G2001.

Definition at line 309 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::model_, bpp::MarkovModulatedSubstitutionModel::updateMatrices_(), and bpp::MarkovModulatedSubstitutionModel::updateRatesModel_().

Referenced by bpp::G2001::fireParameterChanged().

◆ freq()

double bpp::MarkovModulatedSubstitutionModel::freq ( size_t  i) const
inlineoverridevirtualinherited
Returns
Equilibrium frequency associated to character i.
See also
getFrequencies(), getStates()

Implements bpp::TransitionModelInterface.

Definition at line 216 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::freq_.

◆ frequencySet()

const FrequencySetInterface& bpp::MarkovModulatedSubstitutionModel::frequencySet ( ) const
inlineoverridevirtualinherited
Returns
Get the FrequencySet of equilibrium of this model.
Exceptions
Exceptionif no FrequenceSet is associated to this model.

Implements bpp::BranchModelInterface.

Definition at line 238 of file MarkovModulatedSubstitutionModel.h.

◆ generator()

const Matrix<double>& bpp::MarkovModulatedSubstitutionModel::generator ( ) const
inlineoverridevirtualinherited
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. Additionally, the rate_ attribute 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::SubstitutionModelInterface.

Definition at line 201 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::generator_.

◆ getAlphabet()

std::shared_ptr<const Alphabet> bpp::MarkovModulatedSubstitutionModel::getAlphabet ( ) const
inlineoverridevirtualinherited
Returns
A pointer toward the alphabet associated to this model.

Implements bpp::BranchModelInterface.

Definition at line 179 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::model_.

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

◆ getAlphabetStateAsChar()

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

Implements bpp::BranchModelInterface.

Definition at line 189 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::stateMap_.

◆ getAlphabetStateAsInt()

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

Implements bpp::BranchModelInterface.

Definition at line 191 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::stateMap_.

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

◆ getAlphabetStates()

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

Implements bpp::BranchModelInterface.

Definition at line 187 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::stateMap_.

◆ getColumnRightEigenVectors()

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

Implements bpp::SubstitutionModelInterface.

Definition at line 214 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::rightEigenVectors_.

◆ getd2Pij_dt2()

const Matrix< double > & MarkovModulatedSubstitutionModel::getd2Pij_dt2 ( double  t) const
overridevirtualinherited

◆ getdPij_dt()

const Matrix< double > & MarkovModulatedSubstitutionModel::getdPij_dt ( double  t) const
overridevirtualinherited

◆ getEigenValues()

const Vdouble& bpp::MarkovModulatedSubstitutionModel::getEigenValues ( ) const
inlineoverridevirtualinherited
Returns
A vector with all real parts of the eigen values of the generator of this model;

Implements bpp::SubstitutionModelInterface.

Definition at line 207 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::eigenValues_.

◆ getFrequencies()

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

Implements bpp::TransitionModelInterface.

Definition at line 197 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::freq_.

◆ getFrequencies_()

Vdouble& bpp::MarkovModulatedSubstitutionModel::getFrequencies_ ( )
inlineoverrideprotectedvirtualinherited

◆ getIEigenValues()

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

Implements bpp::SubstitutionModelInterface.

Definition at line 208 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::iEigenValues_.

◆ getInitValue()

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

This method is used to initialize likelihoods in recursions. 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::BranchModelInterface.

Definition at line 185 of file MarkovModulatedSubstitutionModel.cpp.

References bpp::MarkovModulatedSubstitutionModel::getAlphabet(), bpp::MarkovModulatedSubstitutionModel::getAlphabetStateAsInt(), bpp::MarkovModulatedSubstitutionModel::model_, bpp::MarkovModulatedSubstitutionModel::nbRates_, and bpp::MarkovModulatedSubstitutionModel::nbStates_.

◆ getModelStates() [1/2]

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

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::BranchModelInterface.

Definition at line 195 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::stateMap_.

◆ getModelStates() [2/2]

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

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::BranchModelInterface.

Definition at line 193 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::stateMap_.

◆ getName()

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

Get the name of the model.

Returns
The name of this model.

Implements bpp::BranchModelInterface.

Definition at line 72 of file TS98.h.

◆ getNumberOfStates()

size_t bpp::MarkovModulatedSubstitutionModel::getNumberOfStates ( ) const
inlineoverridevirtualinherited

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::BranchModelInterface.

Definition at line 181 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::stateMap_.

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

◆ getPij_t()

◆ getRate() [1/2]

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

Get the rate.

Reimplemented from bpp::MarkovModulatedSubstitutionModel.

Definition at line 74 of file TS98.h.

◆ getRate() [2/2]

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 255 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::nbStates_.

◆ getRowLeftEigenVectors()

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

Implements bpp::SubstitutionModelInterface.

Definition at line 213 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::leftEigenVectors_.

◆ getScale()

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

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::SubstitutionModelInterface.

Definition at line 282 of file MarkovModulatedSubstitutionModel.h.

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

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

◆ getStateMap()

std::shared_ptr<const StateMapInterface> bpp::MarkovModulatedSubstitutionModel::getStateMap ( ) const
inlineoverridevirtualinherited
Returns
A shared_ptr to the mapping of model states with alphabet states.

Implements bpp::BranchModelInterface.

Definition at line 185 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::stateMap_.

◆ isDiagonalizable()

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

Implements bpp::SubstitutionModelInterface.

Definition at line 210 of file MarkovModulatedSubstitutionModel.h.

◆ isNonSingular()

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

Implements bpp::SubstitutionModelInterface.

Definition at line 211 of file MarkovModulatedSubstitutionModel.h.

◆ isScalable()

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

◆ Lik_t() [1/2]

virtual const Eigen::VectorXd& bpp::TransitionModelInterface::Lik_t ( const Eigen::VectorXd &  values,
double  t 
) const
pure virtualinherited

This method is used to compute likelihoods in recursions. It computes the probability of a vector given a start state.

Parameters
valuesAn vector of states on the site.
tthe branch length

Implements bpp::BranchModelInterface.

Implemented in bpp::AbstractLkTransitionModel.

◆ Lik_t() [2/2]

const Eigen::VectorXd& bpp::AbstractLkTransitionModel::Lik_t ( const Eigen::VectorXd &  values,
double  t 
) const
inlineoverridevirtualinherited

This method is used to compute likelihoods in recursions. It computes the probability of a vector given a start state.

Parameters
valuesAn vector of states on the site.
tthe branch length

Implements bpp::TransitionModelInterface.

Definition at line 30 of file AbstractSubstitutionModel.h.

References bpp::TransitionModelInterface::getPij_t(), and bpp::AbstractLkTransitionModel::lik_.

◆ nestedModel()

const ReversibleSubstitutionModelInterface& bpp::MarkovModulatedSubstitutionModel::nestedModel ( ) const
inlineinherited

◆ normalize()

void bpp::MarkovModulatedSubstitutionModel::normalize ( )
inlineoverridevirtualinherited

◆ Pij_t()

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

Implements bpp::TransitionModelInterface.

Definition at line 220 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::getPij_t().

◆ Qij()

double bpp::MarkovModulatedSubstitutionModel::Qij ( size_t  i,
size_t  j 
) const
inlineoverridevirtualinherited

A method for computing all necessary matrices.

Returns
The rate in the generator of change from state i to state j.
See also
getStates();

Implements bpp::SubstitutionModelInterface.

Definition at line 218 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::generator_.

◆ setDiagonal()

void MarkovModulatedSubstitutionModel::setDiagonal ( )
overridevirtualinherited

◆ setFreq()

void bpp::MarkovModulatedSubstitutionModel::setFreq ( std::map< int, double > &  frequencies)
inlineoverridevirtualinherited

Set equilibrium frequencies.

Parameters
frequenciesThe map of the frequencies to use.

Implements bpp::TransitionModelInterface.

Definition at line 232 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::model_, and bpp::MarkovModulatedSubstitutionModel::updateMatrices_().

◆ setFreqFromData()

void bpp::MarkovModulatedSubstitutionModel::setFreqFromData ( const SequenceDataInterface data,
double  pseudoCount = 0 
)
inlineoverridevirtualinherited

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::TransitionModelInterface.

Definition at line 226 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::model_, and bpp::MarkovModulatedSubstitutionModel::updateMatrices_().

◆ setNamespace()

void MarkovModulatedSubstitutionModel::setNamespace ( const std::string &  prefix)
overridevirtualinherited

◆ setRate()

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

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

Parameters
ratemust be positive.

Reimplemented from bpp::MarkovModulatedSubstitutionModel.

Definition at line 76 of file TS98.h.

◆ setScalable()

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

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

Implements bpp::SubstitutionModelInterface.

Definition at line 269 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::model_.

◆ setScale()

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

Multiplies the current generator by the given scale.

Parameters
scalethe scale by which the generator is multiplied.

Implements bpp::SubstitutionModelInterface.

Definition at line 289 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::model_, and bpp::MarkovModulatedSubstitutionModel::updateMatrices_().

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

◆ Sij()

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

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

Implements bpp::SubstitutionModelInterface.

Definition at line 217 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::exchangeability_.

◆ stateMap()

const StateMapInterface& bpp::MarkovModulatedSubstitutionModel::stateMap ( ) const
inlineoverridevirtualinherited
Returns
The mapping of model states with alphabet states.

Implements bpp::BranchModelInterface.

Definition at line 183 of file MarkovModulatedSubstitutionModel.h.

References bpp::MarkovModulatedSubstitutionModel::stateMap_.

◆ updateMatrices_()

void MarkovModulatedSubstitutionModel::updateMatrices_ ( )
protectedvirtualinherited

Definition at line 74 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_, RowMatrix< double >::getNumberOfColumns(), bpp::EigenValue< class >::getRealEigenValues(), bpp::MarkovModulatedSubstitutionModel::getScale(), bpp::EigenValue< class >::getV(), bpp::MarkovModulatedSubstitutionModel::iEigenValues_, bpp::MatrixTools::inv(), bpp::MarkovModulatedSubstitutionModel::isScalable(), bpp::MatrixTools::kroneckerMult(), bpp::VectorTools::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_, RowMatrix< double >::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().

◆ updateRatesModel_()

void bpp::TS98::updateRatesModel_ ( )
inlineoverrideprotectedvirtual

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 81 of file TS98.h.

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

Referenced by TS98().

Member Data Documentation

◆ compFreq_

bool bpp::MarkovModulatedSubstitutionModel::compFreq_
protectedinherited

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

Definition at line 110 of file MarkovModulatedSubstitutionModel.h.

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

◆ d2pijt_

RowMatrix<double> bpp::MarkovModulatedSubstitutionModel::d2pijt_
mutableprotectedinherited

◆ dpijt_

RowMatrix<double> bpp::MarkovModulatedSubstitutionModel::dpijt_
mutableprotectedinherited

◆ eigenDecompose_

bool bpp::MarkovModulatedSubstitutionModel::eigenDecompose_
protectedinherited

Tell if the eigen decomposition should be performed.

Definition at line 103 of file MarkovModulatedSubstitutionModel.h.

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

◆ eigenValues_

◆ exchangeability_

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

◆ freq_

◆ generator_

◆ iEigenValues_

Vdouble bpp::MarkovModulatedSubstitutionModel::iEigenValues_
protectedinherited

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

Definition at line 98 of file MarkovModulatedSubstitutionModel.h.

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

◆ leftEigenVectors_

◆ lik_

Eigen::VectorXd bpp::AbstractLkTransitionModel::lik_
mutableprivateinherited

◆ model_

◆ nbRates_

◆ nbStates_

◆ nestedPrefix_

◆ normalizeRateChanges_

bool bpp::MarkovModulatedSubstitutionModel::normalizeRateChanges_
protectedinherited

◆ pijt_

RowMatrix<double> bpp::MarkovModulatedSubstitutionModel::pijt_
mutableprotectedinherited

◆ rates_

◆ ratesExchangeability_

RowMatrix<double> bpp::MarkovModulatedSubstitutionModel::ratesExchangeability_
protectedinherited

◆ ratesFreq_

Vdouble bpp::MarkovModulatedSubstitutionModel::ratesFreq_
protectedinherited

◆ ratesGenerator_

RowMatrix<double> bpp::MarkovModulatedSubstitutionModel::ratesGenerator_
protectedinherited

◆ rightEigenVectors_

◆ stateMap_


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