bpp-phyl3  3.0.0
bpp::SubstitutionModelInterface Class Referenceabstract

Interface for all substitution models. More...

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

+ Inheritance diagram for bpp::SubstitutionModelInterface:
+ Collaboration diagram for bpp::SubstitutionModelInterface:

Public Member Functions

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

Protected Member Functions

virtual VdoublegetFrequencies_ ()=0
 
virtual ParameterListgetParameters_ ()=0
 

Detailed Description

Interface for all substitution models.

A substitution model is based on a Markov generator $Q$, the size of which depends on the alphabet used (4 for nucleotides, 20 for proteins, etc.). Each SubstitutionModel object hence includes a pointer toward an alphabet, and provides a method to retrieve the alphabet used (getAlphabet() method).

What we want from a substitution model is to compute the probabilities of state j at time t geven state j at time 0 ( $P_{i,j}(t)$). Typically, this is computed using the formula

\[ P(t) = e^{r \times t \times Q}, \]

where $P(t)$ is the matrix with all probabilities $P_{i,j}(t)$, and $ r $ the rate. For some models, the $P_{i,j}(t)$'s can be computed analytically.

For more complex models, we need to use a eigen-decomposition of $Q$:

\[ Q = U^{-1} . D . U, \]

where $D = diag(\lambda_i)$ is a diagonal matrix. Hence

\[ P(t) = e^{r \times t \times Q} = U^{-1} . e^{r \times D \times t} . U, \]

where $e^{r \times D \times t} = diag\left(e^{r \times \lambda_i \times t}\right)$ is a diagonal matrix with all terms equal to exp the terms in $D$ multiplied per $ r \times t $. $U$ is the matrix of left eigen vectors (by row), and $U^{-1}$ is the matrix of right eigen vectors (by column). The values in $D$ are the eigen values of $Q$. All $Q,U,U^{-1}$ and $D$ (its diagonal) may be retrieved from the class (getEigenValues(), getRowRightEigenVectors() and getColumnLeftEigenVectors() functions).

First and second order derivatives of $P(t)$ with respect to $t$ can also be retrieved. These methods may be useful for optimization processes. Derivatives may be computed analytically, or using the general formulas:

\[ \frac{\partial P(t)}{\partial t} = U^{-1} . diag\left(r \times \lambda_i \times e^{r \times \lambda_i \times t}\right) . U \]

and

\[ \frac{\partial^2 P(t)}{\partial t^2} = U^{-1} . diag\left(r^2 \times \lambda_i^2 \times e^{r \times \lambda_i \times t}\right) . U \]

If Q is not symmetric, then the eigenvalue matrix D is block diagonal with the real eigenvalues in 1-by-1 blocks and any complex eigenvalues, a + i*b, in 2-by-2 blocks, [a, b; -b, a]. That is, if the complex eigenvalues look like


          a + ib   .        .    .
          .        a - ib   .    .
          .        .        x    .
          .        .        .    y

then D looks like


          a          b      .    .
         -b          a      .    .
          .          .      x    .
          .          .      .    y

and exp(tD) equals


          exp(ta)cos(tb)   exp(ta)sin(tb)  .        .
         -exp(ta)sin(tb)   exp(ta)cos(tb)  .        .
          .                .               exp(tx)  .
          .                .               .        exp(ty)

If U is singular, it cannot be inverted. In this case exp(tQ) is approximated using Taylor decomposition:

\[ P(t) = Id + tQ + \frac{(tQ)^2}{2!} + ... + \frac{(tQ)^n}{n!} + ... \]

To prevent approximation issues, if \f$ max(tQ) \f$ is too high (currently above 0.5), \f$ t \f$ is divided in an ad hoc way (e.g. by \f$ N \f$), and we compute \f$ P(t) = (P(t/N))^N \f$ with a Taylor decomposition for \f$ P(t/N) \f$.

In this case, derivatives according to \f$ t \f$ are computed analytically too.

Definition at line 401 of file SubstitutionModel.h.

Constructor & Destructor Documentation

◆ SubstitutionModelInterface()

bpp::SubstitutionModelInterface::SubstitutionModelInterface ( )
inline

Definition at line 405 of file SubstitutionModel.h.

◆ ~SubstitutionModelInterface()

virtual bpp::SubstitutionModelInterface::~SubstitutionModelInterface ( )
inlinevirtual

Definition at line 406 of file SubstitutionModel.h.

Member Function Documentation

◆ addRateParameter()

◆ alphabet()

◆ clone()

SubstitutionModelInterface* bpp::SubstitutionModelInterface::clone ( ) const
overridepure virtual

Implements bpp::TransitionModelInterface.

Implemented in bpp::ReversibleSubstitutionModelInterface, bpp::AbstractReversibleProteinSubstitutionModel, bpp::AbstractProteinSubstitutionModel, bpp::ProteinReversibleSubstitutionModelInterface, bpp::ProteinSubstitutionModelInterface, bpp::AbstractReversibleNucleotideSubstitutionModel, bpp::AbstractNucleotideSubstitutionModel, bpp::NucleotideReversibleSubstitutionModelInterface, bpp::NucleotideSubstitutionModelInterface, bpp::MarkovModulatedSubstitutionModel, bpp::CodonReversibleSubstitutionModelInterface, bpp::CodonSubstitutionModelInterface, bpp::AbstractKroneckerCodonSubstitutionModel, bpp::AbstractDFPSubstitutionModel, bpp::AbstractCodonSubstitutionModel, bpp::AbstractReversibleSubstitutionModel, bpp::WordSubstitutionModel, bpp::TS98, bpp::RegisterRatesSubstitutionModel, bpp::RE08Codon, bpp::RE08Protein, bpp::RE08Nucleotide, bpp::RE08, bpp::WAG01, bpp::UserProteinSubstitutionModel, bpp::LLG08_UL3::EmbeddedModel, bpp::LLG08_UL2::EmbeddedModel, bpp::LLG08_EX3::EmbeddedModel, bpp::LLG08_EX2::EmbeddedModel, bpp::LLG08_EHO::EmbeddedModel, bpp::LGL08_CAT::EmbeddedModel, bpp::LG08, bpp::JTT92, bpp::JCprot, bpp::DSO78, bpp::Coala, bpp::POMO, bpp::YpR_Gen, bpp::YpR_Sym, bpp::TN93, bpp::T92, bpp::SSR, bpp::RN95s, bpp::RN95, bpp::L95, bpp::K80, bpp::JCnuc, bpp::GTR, bpp::gBGC, bpp::F81, bpp::KroneckerWordSubstitutionModel, bpp::InMixedSubstitutionModel, bpp::G2001, bpp::FromMixtureSubstitutionModel, bpp::EquiprobableSubstitutionModel, bpp::D1WalkSubstitutionModel, bpp::YN98, bpp::TripletSubstitutionModel, bpp::SENCA, bpp::MG94, bpp::KroneckerCodonDistanceSubstitutionModel, bpp::KroneckerCodonDistanceFrequenciesSubstitutionModel, bpp::KCM, bpp::GY94, bpp::DFPDistanceFrequenciesSubstitutionModel, bpp::CodonSameAARateSubstitutionModel, bpp::CodonDistancePhaseFrequenciesSubstitutionModel, bpp::CodonDistanceFrequenciesSubstitutionModel, bpp::CodonAdHocSubstitutionModel, bpp::AnonymousSubstitutionModel, bpp::TwoParameterBinarySubstitutionModel, bpp::LG10_EX_EHO::EmbeddedModel, bpp::HKY85, bpp::F84, bpp::CodonDistanceSubstitutionModel, and bpp::BinarySubstitutionModel.

◆ computeFrequencies() [1/2]

◆ computeFrequencies() [2/2]

virtual void bpp::TransitionModelInterface::computeFrequencies ( bool  yn)
pure virtualinherited

◆ d2Lik_dt2()

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

◆ d2Pij_dt2()

◆ dLik_dt()

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

◆ dPij_dt()

◆ enableEigenDecomposition() [1/2]

virtual bool bpp::SubstitutionModelInterface::enableEigenDecomposition ( )
pure virtual

◆ enableEigenDecomposition() [2/2]

virtual void bpp::SubstitutionModelInterface::enableEigenDecomposition ( bool  yn)
pure virtual

◆ exchangeabilityMatrix()

virtual const Matrix<double>& bpp::SubstitutionModelInterface::exchangeabilityMatrix ( ) const
pure virtual
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$.

Implemented in bpp::MarkovModulatedSubstitutionModel, bpp::InMixedSubstitutionModel, bpp::AbstractTotallyWrappedSubstitutionModel, and bpp::AbstractSubstitutionModel.

Referenced by bpp::Coala::Coala(), bpp::AbstractTotallyWrappedSubstitutionModel::exchangeabilityMatrix(), and bpp::InMixedSubstitutionModel::exchangeabilityMatrix().

◆ freq()

◆ frequencySet()

◆ generator()

virtual const Matrix<double>& bpp::SubstitutionModelInterface::generator ( ) const
pure virtual
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()

Implemented in bpp::MarkovModulatedSubstitutionModel, bpp::InMixedSubstitutionModel, bpp::AbstractTotallyWrappedSubstitutionModel, and bpp::AbstractSubstitutionModel.

Referenced by bpp::OneChangeRegisterTransitionModel::d2Pij_dt2(), bpp::OneChangeRegisterTransitionModel::dPij_dt(), bpp::AbstractTotallyWrappedSubstitutionModel::generator(), bpp::InMixedSubstitutionModel::generator(), bpp::OneChangeRegisterTransitionModel::getd2Pij_dt2(), bpp::OneChangeRegisterTransitionModel::getdPij_dt(), bpp::OneChangeRegisterTransitionModel::getPij_t(), bpp::OneChangeRegisterTransitionModel::Pij_t(), bpp::OneChangeRegisterTransitionModel::updateMatrices_(), and bpp::RegisterRatesSubstitutionModel::updateMatrices_().

◆ getAlphabet()

◆ getAlphabetStateAsChar()

virtual std::string bpp::BranchModelInterface::getAlphabetStateAsChar ( size_t  index) const
pure virtualinherited

◆ getAlphabetStateAsInt()

virtual int bpp::BranchModelInterface::getAlphabetStateAsInt ( size_t  index) const
pure virtualinherited

◆ getAlphabetStates()

virtual const std::vector<int>& bpp::BranchModelInterface::getAlphabetStates ( ) const
pure virtualinherited

◆ getColumnRightEigenVectors()

virtual const Matrix<double>& bpp::SubstitutionModelInterface::getColumnRightEigenVectors ( ) const
pure virtual

◆ getd2Pij_dt2()

◆ getdPij_dt()

◆ getEigenValues()

virtual const Vdouble& bpp::SubstitutionModelInterface::getEigenValues ( ) const
pure virtual

◆ getFrequencies()

◆ getFrequencies_()

◆ getIEigenValues()

virtual const Vdouble& bpp::SubstitutionModelInterface::getIEigenValues ( ) const
pure virtual

◆ getInitValue()

virtual double bpp::BranchModelInterface::getInitValue ( size_t  i,
int  state 
) const
pure virtualinherited

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();

Implemented in bpp::TransitionFromTransitionModel, bpp::RE08, bpp::OneChangeTransitionModel, bpp::OneChangeRegisterTransitionModel, bpp::MultinomialFromTransitionModel, bpp::MarkovModulatedSubstitutionModel, bpp::InMixedSubstitutionModel, bpp::AbstractTotallyWrappedTransitionModel, and bpp::AbstractTransitionModel.

Referenced by bpp::AbstractTotallyWrappedTransitionModel::getInitValue(), bpp::InMixedSubstitutionModel::getInitValue(), bpp::MultinomialFromTransitionModel::getInitValue(), bpp::OneChangeRegisterTransitionModel::getInitValue(), bpp::OneChangeTransitionModel::getInitValue(), and bpp::TransitionFromTransitionModel::getInitValue().

◆ getModelStates() [1/2]

virtual std::vector<size_t> bpp::BranchModelInterface::getModelStates ( const std::string &  code) const
pure virtualinherited

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.

Implemented in bpp::RegisterRatesSubstitutionModel, bpp::MarkovModulatedSubstitutionModel, bpp::AbstractWrappedModel, and bpp::AbstractTransitionModel.

◆ getModelStates() [2/2]

virtual std::vector<size_t> bpp::BranchModelInterface::getModelStates ( int  code) const
pure virtualinherited

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.

Implemented in bpp::RegisterRatesSubstitutionModel, bpp::MarkovModulatedSubstitutionModel, bpp::AbstractWrappedModel, and bpp::AbstractTransitionModel.

Referenced by bpp::SubstitutionModelSet::getModelStates(), and bpp::AbstractWrappedModel::getModelStates().

◆ getName()

virtual std::string bpp::BranchModelInterface::getName ( ) const
pure virtualinherited

Get the name of the model.

Returns
The name of this model.

Implemented in bpp::WordSubstitutionModel, bpp::TS98, bpp::TransitionFromTransitionModel, bpp::RegisterRatesSubstitutionModel, bpp::RE08, bpp::WAG01, bpp::UserProteinSubstitutionModel, bpp::LLG08_UL3, bpp::LLG08_UL3::EmbeddedModel, bpp::LLG08_UL2, bpp::LLG08_UL2::EmbeddedModel, bpp::LLG08_EX3, bpp::LLG08_EX3::EmbeddedModel, bpp::LLG08_EX2, bpp::LLG08_EX2::EmbeddedModel, bpp::LLG08_EHO, bpp::LLG08_EHO::EmbeddedModel, bpp::LGL08_CAT, bpp::LGL08_CAT::EmbeddedModel, bpp::LG08, bpp::JTT92, bpp::JCprot, bpp::DSO78, bpp::Coala, bpp::POMO, bpp::OneChangeTransitionModel, bpp::OneChangeRegisterTransitionModel, bpp::YpR_Gen, bpp::YpR_Sym, bpp::TN93, bpp::T92, bpp::SSR, bpp::RN95s, bpp::RN95, bpp::L95, bpp::K80, bpp::JCnuc, bpp::GTR, bpp::gBGC, bpp::F81, bpp::MultinomialFromTransitionModel, bpp::MixtureOfTransitionModels, bpp::MixtureOfATransitionModel, bpp::KroneckerWordSubstitutionModel, bpp::InMixedSubstitutionModel, bpp::G2001, bpp::FromMixtureSubstitutionModel, bpp::EquiprobableSubstitutionModel, bpp::D1WalkSubstitutionModel, bpp::YNGP_M9, bpp::YNGP_M8, bpp::YNGP_M7, bpp::YNGP_M3, bpp::YNGP_M2, bpp::YNGP_M10, bpp::YNGP_M1, bpp::YN98, bpp::TripletSubstitutionModel, bpp::SENCA, bpp::RELAX, bpp::MG94, bpp::KroneckerCodonDistanceSubstitutionModel, bpp::KroneckerCodonDistanceFrequenciesSubstitutionModel, bpp::KCM, bpp::GY94, bpp::DFPDistanceFrequenciesSubstitutionModel, bpp::DFP07, bpp::CodonSameAARateSubstitutionModel, bpp::CodonDistancePhaseFrequenciesSubstitutionModel, bpp::CodonDistanceFrequenciesSubstitutionModel, bpp::CodonAdHocSubstitutionModel, bpp::AnonymousSubstitutionModel, bpp::AbstractWrappedModel, bpp::TwoParameterBinarySubstitutionModel, bpp::LG10_EX_EHO, bpp::LG10_EX_EHO::EmbeddedModel, bpp::HKY85, bpp::F84, bpp::CodonDistanceSubstitutionModel, and bpp::BinarySubstitutionModel.

Referenced by bpp::YpR::checkModel(), bpp::TransitionMatrixFromModel::compute(), bpp::AbstractWrappedModel::getName(), bpp::FromMixtureSubstitutionModel::getName(), bpp::AbstractSubstitutionModel::getPij_t(), bpp::BppOSubstitutionModelFormat::initialize_(), bpp::MixtureOfSubstitutionModels::MixtureOfSubstitutionModels(), bpp::MixtureOfATransitionModel::model(), bpp::LegacyPhylogeneticsApplicationTools::setSubstitutionModelParametersInitialValuesWithAliases(), bpp::ModelPath::toString(), and bpp::BppOSubstitutionModelFormat::write().

◆ getNumberOfStates()

◆ getPij_t()

◆ getRate()

◆ getRowLeftEigenVectors()

virtual const Matrix<double>& bpp::SubstitutionModelInterface::getRowLeftEigenVectors ( ) const
pure virtual

◆ getScale()

virtual double bpp::SubstitutionModelInterface::getScale ( ) const
pure virtual

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.

Implemented in bpp::MarkovModulatedSubstitutionModel, bpp::InMixedSubstitutionModel, bpp::AbstractTotallyWrappedSubstitutionModel, and bpp::AbstractSubstitutionModel.

Referenced by bpp::AbstractTotallyWrappedSubstitutionModel::getScale(), and bpp::InMixedSubstitutionModel::getScale().

◆ getStateMap()

virtual std::shared_ptr<const StateMapInterface> bpp::BranchModelInterface::getStateMap ( ) const
pure virtualinherited

◆ isDiagonalizable()

virtual bool bpp::SubstitutionModelInterface::isDiagonalizable ( ) const
pure virtual

◆ isNonSingular()

virtual bool bpp::SubstitutionModelInterface::isNonSingular ( ) const
pure virtual

◆ isScalable()

◆ Lik_t()

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.

◆ normalize()

◆ Pij_t()

◆ Qij()

◆ setDiagonal()

virtual void bpp::SubstitutionModelInterface::setDiagonal ( )
pure virtual

◆ setFreq()

virtual void bpp::TransitionModelInterface::setFreq ( std::map< int, double > &  frequencies)
pure virtualinherited

◆ setFreqFromData()

virtual void bpp::TransitionModelInterface::setFreqFromData ( const SequenceDataInterface data,
double  pseudoCount = 0 
)
pure virtualinherited

◆ setRate()

◆ setScalable()

virtual void bpp::SubstitutionModelInterface::setScalable ( bool  scalable)
pure virtual

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

Implemented in bpp::MarkovModulatedSubstitutionModel, bpp::InMixedSubstitutionModel, bpp::AbstractTotallyWrappedSubstitutionModel, and bpp::AbstractSubstitutionModel.

Referenced by bpp::AbstractTotallyWrappedSubstitutionModel::setScalable(), and bpp::InMixedSubstitutionModel::setScalable().

◆ setScale()

virtual void bpp::SubstitutionModelInterface::setScale ( double  scale)
pure virtual

Multiplies the current generator by the given scale.

Parameters
scalethe scale by which the generator is multiplied.

Implemented in bpp::MarkovModulatedSubstitutionModel, bpp::InMixedSubstitutionModel, bpp::AbstractTotallyWrappedSubstitutionModel, and bpp::AbstractSubstitutionModel.

Referenced by bpp::AbstractTotallyWrappedSubstitutionModel::setScale(), and bpp::InMixedSubstitutionModel::setScale().

◆ Sij()

virtual double bpp::SubstitutionModelInterface::Sij ( size_t  i,
size_t  j 
) const
pure virtual

◆ stateMap()

virtual const StateMapInterface& bpp::BranchModelInterface::stateMap ( ) const
pure virtualinherited

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