bpp-phyl  2.2.0
bpp::SubstitutionModel Class Referenceabstract

Interface for all substitution models. More...

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

+ Inheritance diagram for bpp::SubstitutionModel:
+ Collaboration diagram for bpp::SubstitutionModel:

Public Member Functions

 SubstitutionModel ()
 
virtual ~SubstitutionModel ()
 
SubstitutionModelclone () const =0
 
virtual std::string getName () const =0
 Get the name of the model. More...
 
virtual const std::vector< int > & getAlphabetStates () const =0
 
virtual const StateMapgetStateMap () 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 double freq (size_t i) const =0
 
virtual double Qij (size_t i, size_t j) 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 Vdouble & getFrequencies () const =0
 
virtual const Matrix< double > & getGenerator () const =0
 
virtual const Matrix< double > & getExchangeabilityMatrix () const =0
 
virtual double Sij (size_t i, size_t j) const =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 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 Vdouble & getEigenValues () const =0
 
virtual const Vdouble & getIEigenValues () 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 const Alphabet * getAlphabet () const =0
 
virtual size_t getNumberOfStates () const =0
 Get the number of states. More...
 
virtual double getInitValue (size_t i, int state) const =0 throw (IndexOutOfBoundsException, BadIntException)
 
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 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 void setFreqFromData (const SequenceContainer &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)
 Set equilibrium frequencies. More...
 
virtual const FrequenciesSetgetFrequenciesSet () const
 If the model owns a FrequenciesSet, returns a pointer to it, otherwise return 0. More...
 

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 199 of file SubstitutionModel.h.

Constructor & Destructor Documentation

◆ SubstitutionModel()

bpp::SubstitutionModel::SubstitutionModel ( )
inline

Definition at line 203 of file SubstitutionModel.h.

◆ ~SubstitutionModel()

virtual bpp::SubstitutionModel::~SubstitutionModel ( )
inlinevirtual

Definition at line 204 of file SubstitutionModel.h.

Member Function Documentation

◆ addRateParameter()

virtual void bpp::SubstitutionModel::addRateParameter ( )
pure virtual

◆ clone()

SubstitutionModel* bpp::SubstitutionModel::clone ( ) const
pure virtual

Implemented in bpp::ReversibleSubstitutionModel, bpp::AbstractReversibleSubstitutionModel, bpp::YpR_Gen, bpp::AbstractSubstitutionModel, bpp::YpR_Sym, bpp::HKY85, bpp::F84, bpp::JCprot, bpp::MarkovModulatedSubstitutionModel, bpp::T92, bpp::K80, bpp::GTR, bpp::RN95, bpp::MixtureOfSubstitutionModels, bpp::TN93, bpp::RE08, bpp::AbstractCodonSubstitutionModel, bpp::JCnuc, bpp::CodonDistanceFrequenciesSubstitutionModel, bpp::RN95s, bpp::BinarySubstitutionModel, bpp::UserProteinSubstitutionModel, bpp::CodonDistancePhaseFrequenciesSubstitutionModel, bpp::MixtureOfASubstitutionModel, bpp::WAG01, bpp::SSR, bpp::CodonDistanceFitnessPhaseFrequenciesSubstitutionModel, bpp::L95, bpp::DSO78, bpp::JTT92, bpp::LLG08_EX3, bpp::LG08, bpp::LLG08_EHO, bpp::LLG08_UL3, bpp::CodonRateFrequenciesSubstitutionModel, bpp::CodonDistanceSubstitutionModel, bpp::G2001, bpp::LG10_EX_EHO, bpp::LGL08_CAT, bpp::LLG08_EX2, bpp::LLG08_UL2, bpp::AbstractMixedSubstitutionModel, bpp::gBGC, bpp::GY94, bpp::YN98, bpp::TS98, bpp::WordSubstitutionModel, bpp::TripletSubstitutionModel, bpp::MG94, bpp::YNGKP_M1, bpp::CodonRateSubstitutionModel, bpp::YNGKP_M8, bpp::Coala, bpp::YNGKP_M7, bpp::YNGKP_M3, bpp::LLG08_EX3::EmbeddedModel, bpp::YNGKP_M2, bpp::LLG08_EHO::EmbeddedModel, bpp::LLG08_UL3::EmbeddedModel, bpp::LG10_EX_EHO::EmbeddedModel, bpp::LGL08_CAT::EmbeddedModel, bpp::LLG08_EX2::EmbeddedModel, bpp::LLG08_UL2::EmbeddedModel, bpp::AbstractBiblioSubstitutionModel, bpp::AbstractBiblioMixedSubstitutionModel, bpp::CodonSubstitutionModel, bpp::MixedSubstitutionModel, bpp::NucleotideSubstitutionModel, and bpp::ProteinSubstitutionModel.

Referenced by bpp::SubstitutionModelSet::addModel(), bpp::NonHomogeneousSequenceSimulator::NonHomogeneousSequenceSimulator(), and bpp::YpR::operator=().

◆ d2Pij_dt2()

virtual double bpp::SubstitutionModel::d2Pij_dt2 ( size_t  i,
size_t  j,
double  t 
) const
pure virtual
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()

Implemented in bpp::AbstractSubstitutionModel, bpp::MarkovModulatedSubstitutionModel, bpp::F84, bpp::HKY85, bpp::JCprot, bpp::T92, bpp::K80, bpp::RN95, bpp::TN93, bpp::RE08, bpp::JCnuc, bpp::RN95s, bpp::BinarySubstitutionModel, and bpp::AbstractBiblioSubstitutionModel.

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

◆ dPij_dt()

virtual double bpp::SubstitutionModel::dPij_dt ( size_t  i,
size_t  j,
double  t 
) const
pure virtual
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()

Implemented in bpp::AbstractSubstitutionModel, bpp::MarkovModulatedSubstitutionModel, bpp::F84, bpp::HKY85, bpp::JCprot, bpp::T92, bpp::K80, bpp::RN95, bpp::TN93, bpp::RE08, bpp::JCnuc, bpp::RN95s, bpp::BinarySubstitutionModel, and bpp::AbstractBiblioSubstitutionModel.

Referenced by bpp::AbstractBiblioSubstitutionModel::dPij_dt().

◆ enableEigenDecomposition() [1/2]

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

◆ enableEigenDecomposition() [2/2]

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

Tell if eigenValues and Vectors must be computed.

Implemented in bpp::AbstractSubstitutionModel, bpp::MarkovModulatedSubstitutionModel, and bpp::AbstractBiblioSubstitutionModel.

◆ freq()

virtual double bpp::SubstitutionModel::freq ( size_t  i) const
pure virtual

◆ getAlphabet()

virtual const Alphabet* bpp::SubstitutionModel::getAlphabet ( ) const
pure virtual
Returns
Get the alphabet associated to this model.

Implemented in bpp::AbstractSubstitutionModel, bpp::MarkovModulatedSubstitutionModel, and bpp::AbstractBiblioSubstitutionModel.

Referenced by bpp::AAExteriorSubstitutionRegister::AAExteriorSubstitutionRegister(), bpp::AAInteriorSubstitutionRegister::AAInteriorSubstitutionRegister(), bpp::SubstitutionModelSet::addModel(), bpp::DecompositionReward::alphabetIndexHasChanged(), bpp::ComprehensiveSubstitutionRegister::ComprehensiveSubstitutionRegister(), bpp::DecompositionReward::DecompositionReward(), bpp::DecompositionSubstitutionCount::DecompositionSubstitutionCount(), bpp::AbstractMutationProcess::detailedEvolve(), bpp::LaplaceSubstitutionCount::getAllNumbersOfSubstitutions(), bpp::AbstractSubstitutionRegister::getAlphabet(), bpp::AbstractBiblioSubstitutionModel::getAlphabet(), bpp::MarkovModulatedSubstitutionModel::getAlphabet(), bpp::SubstitutionMappingTools::getNormalizationsPerBranch(), bpp::DnDsSubstitutionRegister::getType(), bpp::SelectedSubstitutionRegister::SelectedSubstitutionRegister(), bpp::CategorySubstitutionRegister::setCategories(), bpp::LaplaceSubstitutionCount::setSubstitutionModel(), bpp::UniformizationSubstitutionCount::setSubstitutionModel(), bpp::DecompositionSubstitutionCount::setSubstitutionModel(), bpp::DecompositionReward::setSubstitutionModel(), bpp::UniformizationSubstitutionCount::substitutionRegisterHasChanged(), bpp::DecompositionSubstitutionCount::substitutionRegisterHasChanged(), bpp::UniformizationSubstitutionCount::UniformizationSubstitutionCount(), bpp::YpR::updateMatrices(), and bpp::BppOSubstitutionModelFormat::writeMixed_().

◆ getAlphabetStateAsChar()

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

◆ getAlphabetStateAsInt()

◆ getAlphabetStates()

◆ getColumnRightEigenVectors()

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

◆ getd2Pij_dt2()

◆ getdPij_dt()

◆ getEigenValues()

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

◆ getExchangeabilityMatrix()

virtual const Matrix<double>& bpp::SubstitutionModel::getExchangeabilityMatrix ( ) 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::AbstractSubstitutionModel, bpp::MarkovModulatedSubstitutionModel, and bpp::AbstractBiblioSubstitutionModel.

Referenced by bpp::Coala::Coala(), bpp::AbstractBiblioSubstitutionModel::getExchangeabilityMatrix(), and bpp::MarkovModulatedSubstitutionModel::updateMatrices().

◆ getFrequencies()

◆ getFrequenciesSet()

◆ getGenerator()

virtual const Matrix<double>& bpp::SubstitutionModel::getGenerator ( ) 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. 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()

Implemented in bpp::AbstractSubstitutionModel, bpp::MarkovModulatedSubstitutionModel, and bpp::AbstractBiblioSubstitutionModel.

Referenced by bpp::LaplaceSubstitutionCount::computeCounts(), bpp::UniformizationSubstitutionCount::computeCounts_(), bpp::AbstractBiblioSubstitutionModel::getGenerator(), bpp::SimpleMutationProcess::SimpleMutationProcess(), and bpp::MarkovModulatedSubstitutionModel::updateMatrices().

◆ getIEigenValues()

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

Implemented in bpp::AbstractSubstitutionModel, bpp::MarkovModulatedSubstitutionModel, and bpp::AbstractBiblioSubstitutionModel.

Referenced by bpp::AbstractBiblioSubstitutionModel::getIEigenValues().

◆ getInitValue()

virtual double bpp::SubstitutionModel::getInitValue ( size_t  i,
int  state 
) const
throw (IndexOutOfBoundsException,
BadIntException
)
pure virtual

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

Implemented in bpp::AbstractSubstitutionModel, bpp::MarkovModulatedSubstitutionModel, bpp::RE08, and bpp::AbstractBiblioSubstitutionModel.

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

◆ getModelStates() [1/2]

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

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::AbstractSubstitutionModel, bpp::MarkovModulatedSubstitutionModel, and bpp::AbstractBiblioSubstitutionModel.

Referenced by bpp::AbstractBiblioSubstitutionModel::getModelStates(), bpp::SubstitutionModelSet::getModelStates(), and bpp::MarginalAncestralStateReconstruction::recursiveMarginalAncestralStates().

◆ getModelStates() [2/2]

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

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::AbstractSubstitutionModel, bpp::MarkovModulatedSubstitutionModel, and bpp::AbstractBiblioSubstitutionModel.

◆ getName()

virtual std::string bpp::SubstitutionModel::getName ( ) const
pure virtual

Get the name of the model.

Returns
The name of this model.

Implemented in bpp::YpR_Gen, bpp::YpR_Sym, bpp::F84, bpp::HKY85, bpp::JCprot, bpp::T92, bpp::K80, bpp::RN95, bpp::GTR, bpp::TN93, bpp::RE08, bpp::MixtureOfSubstitutionModels, bpp::JCnuc, bpp::RN95s, bpp::CodonDistanceFrequenciesSubstitutionModel, bpp::BinarySubstitutionModel, bpp::CodonDistancePhaseFrequenciesSubstitutionModel, bpp::UserProteinSubstitutionModel, bpp::CodonDistanceFitnessPhaseFrequenciesSubstitutionModel, bpp::MixtureOfASubstitutionModel, bpp::LLG08_EX3, bpp::WordSubstitutionModel, bpp::LLG08_EHO, bpp::LLG08_UL3, bpp::WAG01, bpp::SSR, bpp::LG10_EX_EHO, bpp::CodonDistanceSubstitutionModel, bpp::L95, bpp::LGL08_CAT, bpp::LLG08_EX2, bpp::LLG08_UL2, bpp::DSO78, bpp::JTT92, bpp::CodonRateFrequenciesSubstitutionModel, bpp::LG08, bpp::YNGKP_M1, bpp::G2001, bpp::gBGC, bpp::YNGKP_M8, bpp::GY94, bpp::YNGKP_M7, bpp::YN98, bpp::YNGKP_M3, bpp::TS98, bpp::TripletSubstitutionModel, bpp::YNGKP_M2, bpp::MG94, bpp::CodonRateSubstitutionModel, bpp::Coala, bpp::LLG08_EX3::EmbeddedModel, bpp::LLG08_EHO::EmbeddedModel, bpp::LLG08_UL3::EmbeddedModel, bpp::LG10_EX_EHO::EmbeddedModel, bpp::LGL08_CAT::EmbeddedModel, bpp::LLG08_EX2::EmbeddedModel, and bpp::LLG08_UL2::EmbeddedModel.

Referenced by bpp::AbstractCodonFitnessSubstitutionModel::AbstractCodonFitnessSubstitutionModel(), bpp::AbstractSubstitutionModel::fireParameterChanged(), bpp::MixtureOfSubstitutionModels::getSubmodelNumbers(), bpp::BppOSubstitutionModelFormat::readMixed_(), bpp::AbstractBiblioSubstitutionModel::updateMatrices(), and bpp::BppOSubstitutionModelFormat::write().

◆ getNumberOfStates()

◆ getPij_t()

◆ getRate()

◆ getRowLeftEigenVectors()

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

◆ getScale()

virtual double bpp::SubstitutionModel::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::AbstractSubstitutionModel, bpp::MarkovModulatedSubstitutionModel, and bpp::AbstractBiblioSubstitutionModel.

Referenced by bpp::AbstractBiblioSubstitutionModel::getScale().

◆ getStateMap()

virtual const StateMap& bpp::SubstitutionModel::getStateMap ( ) const
pure virtual

◆ isDiagonalizable()

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

◆ isNonSingular()

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

◆ Pij_t()

◆ Qij()

◆ setFreq()

◆ setFreqFromData()

virtual void bpp::SubstitutionModel::setFreqFromData ( const SequenceContainer &  data,
double  pseudoCount = 0 
)
pure virtual

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)} \]

Implemented in bpp::AbstractSubstitutionModel, bpp::MarkovModulatedSubstitutionModel, bpp::JCprot, bpp::K80, bpp::RE08, bpp::JCnuc, bpp::UserProteinSubstitutionModel, bpp::AbstractBiblioSubstitutionModel, bpp::WAG01, bpp::DSO78, bpp::JTT92, bpp::LG08, and bpp::Coala.

Referenced by bpp::BppOSubstitutionModelFormat::read(), bpp::AbstractBiblioSubstitutionModel::setFreqFromData(), and bpp::MarkovModulatedSubstitutionModel::setFreqFromData().

◆ setRate()

virtual void bpp::SubstitutionModel::setRate ( double  rate)
pure virtual

◆ setScale()

virtual void bpp::SubstitutionModel::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::AbstractSubstitutionModel, bpp::MarkovModulatedSubstitutionModel, and bpp::AbstractBiblioSubstitutionModel.

Referenced by bpp::AbstractBiblioSubstitutionModel::setScale(), and bpp::MarkovModulatedSubstitutionModel::setScale().

◆ Sij()

virtual double bpp::SubstitutionModel::Sij ( size_t  i,
size_t  j 
) const
pure virtual
Returns
The exchangeability between state i and state j.

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

Implemented in bpp::AbstractSubstitutionModel, bpp::MarkovModulatedSubstitutionModel, and bpp::AbstractBiblioSubstitutionModel.

Referenced by bpp::AbstractBiblioSubstitutionModel::Sij().


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