bpp-phyl  2.4.0
bpp::MixtureOfSubstitutionModels Class Referenceabstract

Substitution models defined as a mixture of several substitution models. More...

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

+ Inheritance diagram for bpp::MixtureOfSubstitutionModels:
+ Collaboration diagram for bpp::MixtureOfSubstitutionModels:

Public Member Functions

 MixtureOfSubstitutionModels (const Alphabet *alpha, std::vector< SubstitutionModel * > vpModel)
 Constructor of a MixtureOfSubstitutionModels, where all the models have rate 1 and equal probability. More...
 
 MixtureOfSubstitutionModels (const Alphabet *alpha, std::vector< SubstitutionModel * > vpModel, Vdouble &vproba, Vdouble &vrate)
 Constructor of a MixtureOfSubstitutionModels. More...
 
 MixtureOfSubstitutionModels (const MixtureOfSubstitutionModels &)
 
MixtureOfSubstitutionModelsoperator= (const MixtureOfSubstitutionModels &)
 
virtual ~MixtureOfSubstitutionModels ()
 
MixtureOfSubstitutionModelsclone () const
 
std::string getName () const
 Get the name of the model. More...
 
const SubstitutionModelgetSubModelWithName (const std::string &name) const
 retrieve a pointer to the submodel with the given name. More...
 
void updateMatrices ()
 Diagonalize the $Q$ matrix, and fill the eigenValues_, iEigenValues_, leftEigenVectors_ and rightEigenVectors_ matrices. More...
 
virtual void setVRates (const Vdouble &vd)
 Sets the rates of the submodels to follow the constraint that the mean rate of the mixture equals rate_. More...
 
Vint getSubmodelNumbers (const std::string &desc) const
 Returns the vector of numbers of the submodels in the mixture that match a description of the parameters numbers. More...
 
void setFreq (std::map< int, double > &)
 applies setFreq to all the models of the mixture and recovers the parameters values. More...
 
virtual size_t getNumberOfModels () const
 returns the number of models in the mixture More...
 
virtual const SubstitutionModelgetNModel (size_t i) const
 Returns a specific model from the mixture. More...
 
virtual SubstitutionModelgetNModel (size_t i)
 
double getNRate (size_t i) const
 Returns the rate of a specific model from the mixture. More...
 
virtual void setRate (double rate)
 Set the rate of the model and the submodels. More...
 
virtual void normalizeVRates ()
 Normalizes the rates of the submodels so that the mean rate of the mixture equals rate_. More...
 
const std::vector< double > & getVRates () const
 Returns the vector of all the rates of the mixture. More...
 
virtual double getNProbability (size_t i) const
 Returns the probability of a specific model from the mixture. More...
 
virtual const std::vector< double > & getProbabilities () const
 Returns the vector of probabilities. More...
 
virtual void setNProbability (size_t i, double prob)
 Sets the probability of a specific model from the mixture. More...
 
double Qij (size_t i, size_t j) const
 This function can not be applied here, so it is defined to prevent wrong usage. More...
 
virtual size_t getNumberOfStates () const
 From SubstitutionModel interface. More...
 
virtual const Matrix< double > & getPij_t (double t) const
 
virtual const Matrix< double > & getdPij_dt (double t) const
 
virtual const Matrix< double > & getd2Pij_dt2 (double t) const
 
virtual const Matrix< double > & getGenerator () const =0
 
const Matrix< double > & getGenerator () const
 
virtual const Matrix< double > & getExchangeabilityMatrix () const =0
 
const Matrix< double > & getExchangeabilityMatrix () const
 
virtual double Sij (size_t i, size_t j) const =0
 
double Sij (size_t i, size_t j) const
 
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...
 
void enableEigenDecomposition (bool yn)
 Set if eigenValues and Vectors must be computed. More...
 
bool enableEigenDecomposition ()
 Tell if eigenValues and Vectors must be computed. More...
 
virtual const VdoublegetEigenValues () const =0
 
const VdoublegetEigenValues () const
 
virtual const VdoublegetIEigenValues () const =0
 
const VdoublegetIEigenValues () const
 
virtual bool isDiagonalizable () const =0
 
bool isDiagonalizable () const
 
virtual bool isNonSingular () const =0
 
bool isNonSingular () const
 
virtual const Matrix< double > & getRowLeftEigenVectors () const =0
 
const Matrix< double > & getRowLeftEigenVectors () const
 
virtual const Matrix< double > & getColumnRightEigenVectors () const =0
 
const Matrix< double > & getColumnRightEigenVectors () const
 
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...
 
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...
 
virtual bool isScalable () const =0
 returns if model is scalable More...
 
virtual bool isScalable () const
 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...
 
double getScale () const
 return scale More...
 
virtual void setScale (double scale)=0
 Multiplies the current generator by the given scale. More...
 
void setScale (double scale)
 Multiplies the current generator by the given scale. More...
 
virtual void normalize ()=0
 Normalize the generator. More...
 
void normalize ()
 normalize the generator More...
 
virtual void setDiagonal ()=0
 set the diagonal of the generator such that sum on each line equals 0. More...
 
void setDiagonal ()
 set the diagonal of the generator such that sum on each line equals 0. More...
 
virtual const std::vector< int > & getAlphabetStates () const =0
 
const std::vector< int > & getAlphabetStates () const
 
virtual const StateMapgetStateMap () const =0
 
const StateMapgetStateMap () const
 
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...
 
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...
 
virtual int getAlphabetStateAsInt (size_t index) const =0
 
int getAlphabetStateAsInt (size_t index) const
 
virtual std::string getAlphabetStateAsChar (size_t index) const =0
 
std::string getAlphabetStateAsChar (size_t index) const
 
virtual double freq (size_t i) const =0
 
virtual double freq (size_t i) const
 
virtual double Pij_t (size_t i, size_t j, double t) const =0
 
virtual double Pij_t (size_t i, size_t j, double t) const
 
virtual double dPij_dt (size_t i, size_t j, double t) const =0
 
virtual double dPij_dt (size_t i, size_t j, double t) const
 
virtual double d2Pij_dt2 (size_t i, size_t j, double t) const =0
 
virtual double d2Pij_dt2 (size_t i, size_t j, double t) const
 
virtual const VdoublegetFrequencies () const =0
 
const VdoublegetFrequencies () const
 
virtual bool computeFrequencies () const =0
 
virtual void computeFrequencies (bool yn)=0
 
bool computeFrequencies () const
 
void computeFrequencies (bool yn)
 
virtual const AlphabetgetAlphabet () const =0
 
const AlphabetgetAlphabet () const
 
virtual double getInitValue (size_t i, int state) const =0
 
double getInitValue (size_t i, int state) const
 
virtual double getRate () const =0
 Get the rate. More...
 
virtual double getRate () const
 The rate of the substitution process. More...
 
virtual void addRateParameter ()=0
 
void addRateParameter ()
 add a "rate" parameter to the model, that handles the overall rate of the process. More...
 
virtual void setFreqFromData (const SequenceContainer &data, double pseudoCount=0)=0
 Set equilibrium frequencies equal to the frequencies estimated from the data. More...
 
void setFreqFromData (const SequenceContainer &data, double pseudoCount=0)
 Set equilibrium frequencies equal to the frequencies estimated from the data. More...
 
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 void setNamespace (const std::string &prefix)=0
 
void setNamespace (const std::string &prefix)
 
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
 
virtual void fireParameterChanged (const ParameterList &parameters)
 Tells the model that a parameter value has changed. More...
 
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

virtual VdoublegetFrequencies_ ()=0
 
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

std::vector< SubstitutionModel * > modelsContainer_
 vector of pointers to SubstitutionModels. More...
 
std::vector< double > vProbas_
 vector of the probabilities of the models More...
 
std::vector< double > vRates_
 vector of the rates of the models. More...
 
const Alphabetalphabet_
 The alphabet relevant to this model. More...
 
std::shared_ptr< const StateMapstateMap_
 The map of model states with alphabet states. More...
 
size_t size_
 The size of the generator, i.e. the number of states. More...
 
bool isScalable_
 If the model is scalable (ie generator can be normalized automatically). More...
 
double rate_
 The rate of the model (default: 1). The generator (and all its vectorial components) is independent of the rate, since it should be normalized. More...
 
RowMatrix< double > generator_
 The generator matrix $Q$ of the model. More...
 
Vdouble freq_
 The vector $\pi_e$ of equilibrium frequencies. More...
 
bool computeFreq_
 if the Frequencies must be computed from the generator More...
 
RowMatrix< double > exchangeability_
 The exchangeability matrix $S$ of the model, defined as $ S_{ij}=\frac{Q_{ij}}{\pi_j}$. When the model is reversible, this matrix is symetric. More...
 
RowMatrix< double > pijt_
 These ones are for bookkeeping: More...
 
RowMatrix< double > dpijt_
 
RowMatrix< double > d2pijt_
 
bool eigenDecompose_
 Tell if the eigen decomposition should be performed. More...
 
Vdouble eigenValues_
 The vector of eigen values. More...
 
Vdouble iEigenValues_
 The vector of the imaginary part of the eigen values. More...
 
bool isDiagonalizable_
 boolean value for diagonalizability in R of the generator_ More...
 
RowMatrix< double > rightEigenVectors_
 The $U^-1$ matrix made of right eigen vectors (by column). More...
 
bool isNonSingular_
 boolean value for non-singularity of rightEigenVectors_ More...
 
RowMatrix< double > leftEigenVectors_
 The $U$ matrix made of left eigen vectors (by row) if rightEigenVectors_ is non-singular. More...
 
std::vector< RowMatrix< double > > vPowGen_
 vector of the powers of generator_ for Taylor development (if rightEigenVectors_ is singular). More...
 
RowMatrix< double > tmpMat_
 For computational issues. More...
 

Detailed Description

Substitution models defined as a mixture of several substitution models.

Author
Laurent Gu├ęguen

All the models can be of different types (for example T92 or GY94), and each model has a specific probability and rate.

The probabilities and rates of the models are independent parameters, handled directly, under the constraint that the expectation of the rates on the distribution of the models must equal one.

If there are $n$ models, $p_i$ is the probability of model i ( $\sum_{i=1}^{n} p_i = 1$) and the probabilities are defined by relative probabilities parameters $rp_i$ (called "relprobai") with:

\[ 1 <= i < n, p_i = (1-rp_1)*(1-rp_2)...(1-rp_{i-1})*rp_{i} \]

\[ p_n = (1-rp_1)*(1-rp_2)...(1-rp_{n-1}) \]

and

\[ \forall 1 <= i < n, rp_i = \frac{p_i}{1-(p_1+...+p_{i-1})} \]

where $p_i$ stands for the probability of model $i$.

If there are $n$ models, $\rho_i$ is the rate and $p_i$ is the probability of model i ( $\sum_{i=1}^{n} p_i * \rho_i = 1$), the rates are defined by relative rates parameters $r_i$ (called "relratei") with:

\[ 1 <= i < n, \rho_i = (1-r_1)*(1-r_2)...(1-r_{i-1})*\frac{r_{i}}{p_i} \]

\[ \rho_n = \frac{(1-r_1)*(1-r_2)*...*(1-r_{n-1})}{p_n} \]

and

\[ \forall 1 <= i < n, r_i = \frac{\rho_i*p_i}{1-(p_1*\rho_1+...+p_{i-1}*\rho_{i-1})} < 1. \]

For example:

Mixture(model1=HKY85(kappa=3), model2=T92(theta=0.1), model2=L95(gamma=2), relrate1=0.2, relrate2=0.9, relproba1=0.1,relproba2=0.8)

define a model as a mixture of 3 different models: HKY85 has probability 0.1 and rate 2, T92 has probability 0.4 and rate 1.8, and L95 has probability 0.5 and rate 0.16.

The parameters are named "Mixture.relrate1", "Mixture.relrate2", "Mixture.relproba1", "Mixture.relproba2"... in addition to the parameters of the submodels that are prefixed by "Mixture.i_", where i is the order of the model.

Definition at line 117 of file MixtureOfSubstitutionModels.h.

Constructor & Destructor Documentation

MixtureOfSubstitutionModels::MixtureOfSubstitutionModels ( const Alphabet alpha,
std::vector< SubstitutionModel * >  vpModel 
)

Constructor of a MixtureOfSubstitutionModels, where all the models have rate 1 and equal probability.

Parameters
alphapointer to the Alphabet
vpModelvector of pointers to SubstitutionModels. All the SubstitutionModels are owned by the instance.
Warning
providing a vpModel with size 0 will generate a segmentation fault!

Definition at line 49 of file MixtureOfSubstitutionModels.cpp.

References bpp::AbstractParameterAliasable::addParameter_(), bpp::AbstractParameterAliasable::addParameters_(), bpp::AbstractParameterAliasable::getNamespace(), bpp::AbstractParameterAliasable::getParameters(), bpp::AbstractMixedSubstitutionModel::modelsContainer_, bpp::Parameter::PROP_CONSTRAINT_EX, bpp::TextTools::toString(), updateMatrices(), bpp::AbstractMixedSubstitutionModel::vProbas_, and bpp::AbstractMixedSubstitutionModel::vRates_.

Referenced by clone().

MixtureOfSubstitutionModels::MixtureOfSubstitutionModels ( const Alphabet alpha,
std::vector< SubstitutionModel * >  vpModel,
Vdouble vproba,
Vdouble vrate 
)

Constructor of a MixtureOfSubstitutionModels.

Parameters
alphapointer to the Alphabet
vpModelvector of pointers to SubstitutionModels. All the SubstitutionModels are owned by the instance.
vprobavector of the probabilities of the models
vratevector of the rates of the models
Warning
providing a vpModel with size 0 will generate a segmentation fault!

See above the constraints on the rates and the probabilities of the vectors.

Definition at line 98 of file MixtureOfSubstitutionModels.cpp.

References bpp::AbstractParameterAliasable::addParameter_(), bpp::AbstractParameterAliasable::addParameters_(), bpp::AbstractParameterAliasable::getNamespace(), bpp::AbstractParameterAliasable::getParameters(), bpp::AbstractMixedSubstitutionModel::modelsContainer_, bpp::Parameter::PROP_CONSTRAINT_EX, bpp::NumConstants::SMALL(), bpp::TextTools::toString(), updateMatrices(), bpp::AbstractMixedSubstitutionModel::vProbas_, and bpp::AbstractMixedSubstitutionModel::vRates_.

MixtureOfSubstitutionModels::MixtureOfSubstitutionModels ( const MixtureOfSubstitutionModels msm)

Definition at line 179 of file MixtureOfSubstitutionModels.cpp.

MixtureOfSubstitutionModels::~MixtureOfSubstitutionModels ( )
virtual

Definition at line 192 of file MixtureOfSubstitutionModels.cpp.

Member Function Documentation

void AbstractSubstitutionModel::addRateParameter ( )
virtualinherited
MixtureOfSubstitutionModels* bpp::MixtureOfSubstitutionModels::clone ( ) const
inlinevirtual
virtual void bpp::TransitionModel::computeFrequencies ( bool  yn)
pure virtualinherited
void bpp::AbstractSubstitutionModel::computeFrequencies ( bool  yn)
inlinevirtualinherited
Returns
Set if equilibrium frequencies should be computed from the generator

Implements bpp::TransitionModel.

Definition at line 269 of file AbstractSubstitutionModel.h.

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

Reimplemented in bpp::HKY85, bpp::F84, bpp::JCprot, bpp::T92, bpp::K80, bpp::RE08, bpp::RN95, bpp::TN93, bpp::JCnuc, bpp::RN95s, and bpp::BinarySubstitutionModel.

Definition at line 299 of file AbstractSubstitutionModel.h.

References bpp::AbstractSubstitutionModel::getd2Pij_dt2(), bpp::AbstractSubstitutionModel::getInitValue(), bpp::AbstractSubstitutionModel::setFreq(), and bpp::AbstractSubstitutionModel::setFreqFromData().

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

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

Reimplemented in bpp::HKY85, bpp::F84, bpp::JCprot, bpp::T92, bpp::K80, bpp::RE08, bpp::RN95, bpp::TN93, bpp::JCnuc, bpp::RN95s, and bpp::BinarySubstitutionModel.

Definition at line 298 of file AbstractSubstitutionModel.h.

References bpp::AbstractSubstitutionModel::getdPij_dt().

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

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

Set if eigenValues and Vectors must be computed.

Implements bpp::SubstitutionModel.

Definition at line 307 of file AbstractSubstitutionModel.h.

virtual bool bpp::SubstitutionModel::enableEigenDecomposition ( )
pure virtualinherited
virtual void bpp::AbstractSubstitutionModel::fireParameterChanged ( const ParameterList parameters)
inlinevirtualinherited

Tells the model that a parameter value has changed.

This updates the matrices consequently.

Reimplemented from bpp::AbstractParameterAliasable.

Reimplemented in bpp::RegisterRatesSubstitutionModel, bpp::JCprot, bpp::KroneckerCodonDistanceFrequenciesSubstitutionModel, bpp::YpR, bpp::RE08, bpp::KroneckerCodonDistanceSubstitutionModel, bpp::CodonDistanceFrequenciesSubstitutionModel, bpp::CodonAdHocSubstitutionModel, bpp::CodonDistancePhaseFrequenciesSubstitutionModel, bpp::UserProteinSubstitutionModel, bpp::JTT92, bpp::WAG01, bpp::DSO78, bpp::CodonDistanceSubstitutionModel, bpp::LG08, bpp::gBGC, and bpp::SENCA.

Definition at line 316 of file AbstractSubstitutionModel.h.

References bpp::AbstractSubstitutionModel::addRateParameter(), bpp::AbstractParameterAliasable::fireParameterChanged(), bpp::AbstractParameterAliasable::getNamespace(), bpp::ParameterList::getParameterValue(), bpp::ParameterList::hasParameter(), bpp::ParameterList::size(), and bpp::AbstractSubstitutionModel::updateMatrices().

Referenced by bpp::SENCA::fireParameterChanged(), bpp::gBGC::fireParameterChanged(), bpp::LG08::fireParameterChanged(), bpp::CodonDistanceSubstitutionModel::fireParameterChanged(), bpp::DSO78::fireParameterChanged(), bpp::WAG01::fireParameterChanged(), bpp::JTT92::fireParameterChanged(), bpp::UserProteinSubstitutionModel::fireParameterChanged(), bpp::CodonDistancePhaseFrequenciesSubstitutionModel::fireParameterChanged(), bpp::CodonAdHocSubstitutionModel::fireParameterChanged(), bpp::CodonDistanceFrequenciesSubstitutionModel::fireParameterChanged(), bpp::KroneckerCodonDistanceSubstitutionModel::fireParameterChanged(), bpp::KroneckerCodonDistanceFrequenciesSubstitutionModel::fireParameterChanged(), and bpp::JCprot::fireParameterChanged().

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

Implements bpp::TransitionModel.

Definition at line 293 of file AbstractSubstitutionModel.h.

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

Implemented in bpp::RE08Protein, bpp::AbstractSubstitutionModel, bpp::RE08Nucleotide, bpp::MarkovModulatedSubstitutionModel, bpp::RegisterRatesSubstitutionModel, bpp::AbstractReversibleNucleotideSubstitutionModel, bpp::AbstractReversibleProteinSubstitutionModel, bpp::AbstractNucleotideSubstitutionModel, bpp::AbstractProteinSubstitutionModel, bpp::AbstractWrappedModel, bpp::NucleotideSubstitutionModel, and bpp::ProteinSubstitutionModel.

Referenced by bpp::SubstitutionModelSet::addModel(), bpp::DecompositionReward::alphabetIndexHasChanged(), bpp::YpR::check_model(), bpp::ComprehensiveSubstitutionRegister::ComprehensiveSubstitutionRegister(), bpp::SubstitutionMappingTools::computeCountsPerSitePerBranchPerType(), bpp::SubstitutionMappingTools::computeCountsPerSitePerType(), bpp::SubstitutionMappingTools::computeCountsPerTypePerBranch(), bpp::SubstitutionModelSetTools::createHomogeneousModelSet(), bpp::SubstitutionModelSetTools::createNonHomogeneousModelSet(), bpp::DecompositionReward::DecompositionReward(), bpp::DecompositionSubstitutionCount::DecompositionSubstitutionCount(), bpp::LaplaceSubstitutionCount::getAllNumbersOfSubstitutions(), bpp::AbstractWrappedModel::getAlphabet(), bpp::AbstractSubstitutionRegister::getAlphabet(), bpp::MarkovModulatedSubstitutionModel::getAlphabet(), bpp::MarkovModulatedSubstitutionModel::getInitValue(), bpp::SubstitutionMappingTools::getNormalizationsPerBranch(), bpp::PhylogeneticsApplicationTools::getSubstitutionRegister(), bpp::DRASRTreeLikelihoodData::initLikelihoods(), bpp::DRASDRTreeLikelihoodData::initLikelihoods(), bpp::SENCA::setFreq(), bpp::AbstractHomogeneousTreeLikelihood::setModel(), bpp::DecompositionSubstitutionCount::setSubstitutionModel(), bpp::DecompositionReward::setSubstitutionModel(), bpp::LaplaceSubstitutionCount::setSubstitutionModel(), bpp::UniformizationSubstitutionCount::setSubstitutionModel(), bpp::DecompositionSubstitutionCount::substitutionRegisterHasChanged(), bpp::UniformizationSubstitutionCount::substitutionRegisterHasChanged(), bpp::TwoTreeLikelihood::TwoTreeLikelihood(), bpp::UniformizationSubstitutionCount::UniformizationSubstitutionCount(), bpp::YpR::updateMatrices(), bpp::AAInteriorSubstitutionRegister::updateMatrix_(), bpp::AAExteriorSubstitutionRegister::updateMatrix_(), and bpp::BppOSubstitutionModelFormat::writeMixed_().

virtual std::string bpp::TransitionModel::getAlphabetStateAsChar ( size_t  index) const
pure virtualinherited
std::string bpp::AbstractSubstitutionModel::getAlphabetStateAsChar ( size_t  index) const
inlinevirtualinherited
Parameters
indexThe model state.
Returns
The corresponding alphabet state as character code.

Implements bpp::TransitionModel.

Reimplemented in bpp::RegisterRatesSubstitutionModel.

Definition at line 257 of file AbstractSubstitutionModel.h.

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

Implements bpp::TransitionModel.

Reimplemented in bpp::RegisterRatesSubstitutionModel.

Definition at line 259 of file AbstractSubstitutionModel.h.

Referenced by bpp::AbstractCodonSubstitutionModel::completeMatrices(), and bpp::AbstractSubstitutionModel::getInitValue().

const std::vector<int>& bpp::AbstractSubstitutionModel::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.

Reimplemented in bpp::RegisterRatesSubstitutionModel.

Definition at line 255 of file AbstractSubstitutionModel.h.

const Matrix<double>& bpp::AbstractSubstitutionModel::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 291 of file AbstractSubstitutionModel.h.

References bpp::AbstractSubstitutionModel::rightEigenVectors_.

const Matrix< double > & AbstractMixedSubstitutionModel::getd2Pij_dt2 ( double  t) const
virtualinherited
Returns
All second order derivatives of the probability of change from state i to state j with respect to time t, at time t.
See also
d2Pij_dt2()

Reimplemented from bpp::AbstractSubstitutionModel.

Definition at line 170 of file AbstractMixedSubstitutionModel.cpp.

References bpp::AbstractSubstitutionModel::d2pijt_, bpp::AbstractMixedSubstitutionModel::getNumberOfStates(), bpp::AbstractMixedSubstitutionModel::modelsContainer_, and bpp::AbstractMixedSubstitutionModel::vProbas_.

Referenced by bpp::AbstractMixedSubstitutionModel::Qij().

const Matrix< double > & AbstractMixedSubstitutionModel::getdPij_dt ( double  t) const
virtualinherited
Returns
Get all first order derivatives of the probability of change from state i to state j with respect to time t, at time t.
See also
dPij_dt()

Reimplemented from bpp::AbstractSubstitutionModel.

Definition at line 144 of file AbstractMixedSubstitutionModel.cpp.

References bpp::AbstractSubstitutionModel::dpijt_, bpp::AbstractMixedSubstitutionModel::getNumberOfStates(), bpp::AbstractMixedSubstitutionModel::modelsContainer_, and bpp::AbstractMixedSubstitutionModel::vProbas_.

Referenced by bpp::AbstractMixedSubstitutionModel::Qij().

const Vdouble& bpp::AbstractSubstitutionModel::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 281 of file AbstractSubstitutionModel.h.

References bpp::AbstractSubstitutionModel::eigenValues_.

virtual const Vdouble& bpp::SubstitutionModel::getEigenValues ( ) const
pure virtualinherited
const Matrix<double>& bpp::AbstractSubstitutionModel::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 273 of file AbstractSubstitutionModel.h.

References bpp::AbstractSubstitutionModel::exchangeability_.

virtual const Matrix<double>& bpp::SubstitutionModel::getExchangeabilityMatrix ( ) const
pure virtualinherited
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::AbstractTotallyWrappedSubstitutionModel, bpp::MarkovModulatedSubstitutionModel, and bpp::InMixedSubstitutionModel.

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

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

Implements bpp::TransitionModel.

Definition at line 265 of file AbstractSubstitutionModel.h.

References bpp::AbstractSubstitutionModel::freq_.

Referenced by bpp::SENCA::setFreq().

const Matrix<double>& bpp::AbstractSubstitutionModel::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 271 of file AbstractSubstitutionModel.h.

References bpp::AbstractSubstitutionModel::generator_.

Referenced by bpp::AbstractKroneckerWordSubstitutionModel::fillBasicGenerator().

virtual const Matrix<double>& bpp::SubstitutionModel::getGenerator ( ) const
pure virtualinherited
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::AbstractTotallyWrappedSubstitutionModel, bpp::MarkovModulatedSubstitutionModel, and bpp::InMixedSubstitutionModel.

Referenced by bpp::LaplaceSubstitutionCount::computeCounts(), bpp::UniformizationSubstitutionCount::computeCounts_(), bpp::OneChangeRegisterTransitionModel::d2Pij_dt2(), bpp::OneChangeRegisterTransitionModel::dPij_dt(), bpp::OneChangeRegisterTransitionModel::getd2Pij_dt2(), bpp::OneChangeRegisterTransitionModel::getdPij_dt(), bpp::InMixedSubstitutionModel::getGenerator(), bpp::OneChangeRegisterTransitionModel::getPij_t(), bpp::OneChangeRegisterTransitionModel::Pij_t(), bpp::SimpleMutationProcess::SimpleMutationProcess(), bpp::OneChangeRegisterTransitionModel::updateMatrices(), bpp::RegisterRatesSubstitutionModel::updateMatrices(), and bpp::MarkovModulatedSubstitutionModel::updateMatrices().

const Vdouble& bpp::AbstractSubstitutionModel::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 283 of file AbstractSubstitutionModel.h.

References bpp::AbstractSubstitutionModel::iEigenValues_.

virtual const Vdouble& bpp::SubstitutionModel::getIEigenValues ( ) const
pure virtualinherited
virtual double bpp::TransitionModel::getInitValue ( size_t  i,
int  state 
) const
pure 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();

Implemented in bpp::AbstractSubstitutionModel, bpp::MarkovModulatedSubstitutionModel, bpp::RE08, bpp::OneChangeRegisterTransitionModel, bpp::InMixedSubstitutionModel, bpp::AbstractTotallyWrappedModel, and bpp::OneChangeTransitionModel.

Referenced by bpp::OneChangeTransitionModel::getInitValue(), bpp::AbstractTotallyWrappedModel::getInitValue(), bpp::InMixedSubstitutionModel::getInitValue(), bpp::DRASRTreeLikelihoodData::initLikelihoods(), bpp::DRASDRTreeLikelihoodData::initLikelihoods(), bpp::DRASRTreeLikelihoodData::initLikelihoodsWithPatterns(), and bpp::TwoTreeLikelihood::initTreeLikelihoods().

double AbstractSubstitutionModel::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.

Reimplemented in bpp::RE08.

Definition at line 541 of file AbstractSubstitutionModel.cpp.

References bpp::AbstractSubstitutionModel::alphabet_, bpp::Alphabet::getAlias(), bpp::AbstractSubstitutionModel::getAlphabetStateAsInt(), bpp::Alphabet::intToChar(), and bpp::AbstractSubstitutionModel::size_.

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

virtual std::vector<size_t> bpp::TransitionModel::getModelStates ( int  code) const
pure virtualinherited
virtual std::vector<size_t> bpp::TransitionModel::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::AbstractSubstitutionModel, bpp::MarkovModulatedSubstitutionModel, bpp::RegisterRatesSubstitutionModel, and bpp::AbstractWrappedModel.

std::vector<size_t> bpp::AbstractSubstitutionModel::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.

Reimplemented in bpp::RegisterRatesSubstitutionModel.

Definition at line 261 of file AbstractSubstitutionModel.h.

std::vector<size_t> bpp::AbstractSubstitutionModel::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.

Reimplemented in bpp::RegisterRatesSubstitutionModel.

Definition at line 263 of file AbstractSubstitutionModel.h.

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

Get the name of the model.

Returns
The name of this model.

Implements bpp::TransitionModel.

Definition at line 162 of file MixtureOfSubstitutionModels.h.

References getSubmodelNumbers(), getSubModelWithName(), setFreq(), setVRates(), and updateMatrices().

Referenced by getSubModelWithName().

virtual const SubstitutionModel* bpp::AbstractMixedSubstitutionModel::getNModel ( size_t  i) const
inlinevirtualinherited
virtual SubstitutionModel* bpp::AbstractMixedSubstitutionModel::getNModel ( size_t  i)
inlinevirtualinherited

Implements bpp::MixedSubstitutionModel.

Definition at line 127 of file AbstractMixedSubstitutionModel.h.

virtual double bpp::AbstractMixedSubstitutionModel::getNProbability ( size_t  i) const
inlinevirtualinherited

Returns the probability of a specific model from the mixture.

Implements bpp::MixedSubstitutionModel.

Definition at line 179 of file AbstractMixedSubstitutionModel.h.

double bpp::AbstractMixedSubstitutionModel::getNRate ( size_t  i) const
inlinevirtualinherited
virtual size_t bpp::AbstractMixedSubstitutionModel::getNumberOfModels ( ) const
inlinevirtualinherited
const Matrix< double > & AbstractMixedSubstitutionModel::getPij_t ( double  t) const
virtualinherited
virtual const std::vector<double>& bpp::AbstractMixedSubstitutionModel::getProbabilities ( ) const
inlinevirtualinherited

Returns the vector of probabilities.

Implements bpp::MixedSubstitutionModel.

Definition at line 189 of file AbstractMixedSubstitutionModel.h.

References bpp::AbstractMixedSubstitutionModel::vProbas_.

double AbstractSubstitutionModel::getRate ( ) const
virtualinherited

The rate of the substitution process.

Implements bpp::TransitionModel.

Reimplemented in bpp::RegisterRatesSubstitutionModel.

Definition at line 643 of file AbstractSubstitutionModel.cpp.

References bpp::AbstractSubstitutionModel::rate_.

Referenced by bpp::AbstractSubstitutionModel::isScalable().

const Matrix<double>& bpp::AbstractSubstitutionModel::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 289 of file AbstractSubstitutionModel.h.

References bpp::AbstractSubstitutionModel::leftEigenVectors_.

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

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, bpp::AbstractTotallyWrappedSubstitutionModel, and bpp::InMixedSubstitutionModel.

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

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

Implements bpp::TransitionModel.

Reimplemented in bpp::RegisterRatesSubstitutionModel.

Definition at line 253 of file AbstractSubstitutionModel.h.

References bpp::AbstractSubstitutionModel::stateMap_.

Vint MixtureOfSubstitutionModels::getSubmodelNumbers ( const std::string &  desc) const
virtual

Returns the vector of numbers of the submodels in the mixture that match a description of the parameters numbers.

Parameters
descis the description of the class indexes of the mixed parameters. Syntax is like: kappa_1,gamma_3,delta_2

Implements bpp::MixedSubstitutionModel.

Definition at line 308 of file MixtureOfSubstitutionModels.cpp.

References bpp::TransitionModel::getName(), bpp::AbstractMixedSubstitutionModel::getNModel(), and bpp::AbstractMixedSubstitutionModel::getNumberOfModels().

Referenced by getName().

const SubstitutionModel * MixtureOfSubstitutionModels::getSubModelWithName ( const std::string &  name) const
virtual

retrieve a pointer to the submodel with the given name.

Return Null if not found.

Implements bpp::MixedSubstitutionModel.

Definition at line 195 of file MixtureOfSubstitutionModels.cpp.

References getName(), bpp::AbstractMixedSubstitutionModel::getNModel(), and bpp::AbstractMixedSubstitutionModel::getNumberOfModels().

Referenced by getName().

const std::vector<double>& bpp::AbstractMixedSubstitutionModel::getVRates ( ) const
inlinevirtualinherited

Returns the vector of all the rates of the mixture.

Implements bpp::MixedSubstitutionModel.

Definition at line 170 of file AbstractMixedSubstitutionModel.h.

References bpp::AbstractMixedSubstitutionModel::vRates_.

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

Implements bpp::SubstitutionModel.

Definition at line 285 of file AbstractSubstitutionModel.h.

References bpp::AbstractSubstitutionModel::isDiagonalizable_.

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

Implements bpp::SubstitutionModel.

Definition at line 287 of file AbstractSubstitutionModel.h.

References bpp::AbstractSubstitutionModel::isNonSingular_.

virtual bool bpp::SubstitutionModel::isNonSingular ( ) const
pure virtualinherited
virtual bool bpp::SubstitutionModel::isScalable ( ) const
pure virtualinherited
void AbstractMixedSubstitutionModel::normalizeVRates ( )
virtualinherited
MixtureOfSubstitutionModels & MixtureOfSubstitutionModels::operator= ( const MixtureOfSubstitutionModels msm)
virtual double bpp::AbstractSubstitutionModel::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.

Reimplemented in bpp::HKY85, bpp::F84, bpp::JCprot, bpp::T92, bpp::K80, bpp::RE08, bpp::RN95, bpp::TN93, bpp::JCnuc, bpp::RN95s, and bpp::BinarySubstitutionModel.

Definition at line 297 of file AbstractSubstitutionModel.h.

References bpp::AbstractSubstitutionModel::getPij_t().

Referenced by bpp::JCprot::Pij_t().

double bpp::AbstractMixedSubstitutionModel::Qij ( size_t  i,
size_t  j 
) const
inlinevirtualinherited
virtual void bpp::SubstitutionModel::setDiagonal ( )
pure virtualinherited
void MixtureOfSubstitutionModels::setFreq ( std::map< int, double > &  m)
virtual
void AbstractSubstitutionModel::setFreqFromData ( const SequenceContainer data,
double  pseudoCount = 0 
)
virtualinherited

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.

Reimplemented in bpp::JCprot, bpp::K80, bpp::RE08, bpp::JCnuc, bpp::UserProteinSubstitutionModel, bpp::WAG01, bpp::JTT92, bpp::DSO78, bpp::LG08, and bpp::Coala.

Definition at line 559 of file AbstractSubstitutionModel.cpp.

References bpp::SequenceContainerTools::getCounts(), bpp::AbstractSubstitutionModel::setFreq(), and bpp::AbstractSubstitutionModel::size_.

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

virtual void bpp::AbstractMixedSubstitutionModel::setNProbability ( size_t  i,
double  prob 
)
inlinevirtualinherited

Sets the probability of a specific model from the mixture.

Implements bpp::MixedSubstitutionModel.

Definition at line 197 of file AbstractMixedSubstitutionModel.h.

void AbstractMixedSubstitutionModel::setRate ( double  rate)
virtualinherited
void bpp::AbstractSubstitutionModel::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 376 of file AbstractSubstitutionModel.h.

virtual void bpp::SubstitutionModel::setScalable ( bool  scalable)
pure virtualinherited

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

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

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

virtual void bpp::SubstitutionModel::setScale ( double  scale)
pure virtualinherited

Multiplies the current generator by the given scale.

Parameters
scalethe scale by which the generator is multiplied.

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

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

void MixtureOfSubstitutionModels::setVRates ( const Vdouble vd)
virtual

Sets the rates of the submodels to follow the constraint that the mean rate of the mixture equals rate_.

Parameters
vda vector of positive values such that the rates of the respective submodels are in the same proportions (ie this vector does not need to be normalized).

Reimplemented from bpp::AbstractMixedSubstitutionModel.

Definition at line 289 of file MixtureOfSubstitutionModels.cpp.

References bpp::AbstractMixedSubstitutionModel::modelsContainer_, bpp::AbstractParameterAliasable::setParameterValue(), bpp::AbstractMixedSubstitutionModel::setVRates(), bpp::TextTools::toString(), bpp::AbstractMixedSubstitutionModel::vProbas_, and bpp::AbstractMixedSubstitutionModel::vRates_.

Referenced by getName().

double bpp::AbstractSubstitutionModel::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 275 of file AbstractSubstitutionModel.h.

References bpp::AbstractSubstitutionModel::exchangeability_, bpp::AbstractSubstitutionModel::getd2Pij_dt2(), bpp::AbstractSubstitutionModel::getdPij_dt(), and bpp::AbstractSubstitutionModel::getPij_t().

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

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

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

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

void MixtureOfSubstitutionModels::updateMatrices ( )
virtual

Diagonalize the $Q$ matrix, and fill the eigenValues_, iEigenValues_, leftEigenVectors_ and rightEigenVectors_ matrices.

The generator_ matrix and freq_ vector must be initialized.

Eigen values and vectors are computed from the generator and assigned to the eigenValues_ for the real part, iEigenValues_ for the imaginary part, rightEigenVectors_ and leftEigenVectors_ variables. isDiagonalizable_ checks if the generator_ is diagonalizable in R.

The optional rate parameter is not taken into account in this method to prevent unnecessary computation.

!! Here there is no normalization of the generator.

Now check inversion and diagonalization

Reimplemented from bpp::AbstractSubstitutionModel.

Definition at line 206 of file MixtureOfSubstitutionModels.cpp.

References bpp::AbstractSubstitutionModel::freq_, bpp::AbstractMixedSubstitutionModel::getNumberOfStates(), bpp::AbstractParameterAliasable::getParameters(), bpp::AbstractParameterAliasable::getParameterValue(), bpp::AbstractMixedSubstitutionModel::modelsContainer_, bpp::AbstractSubstitutionModel::rate_, bpp::NumConstants::SMALL(), bpp::TextTools::toString(), bpp::AbstractMixedSubstitutionModel::vProbas_, and bpp::AbstractMixedSubstitutionModel::vRates_.

Referenced by getName(), and MixtureOfSubstitutionModels().

Member Data Documentation

bool bpp::AbstractSubstitutionModel::computeFreq_
protectedinherited
bool bpp::AbstractSubstitutionModel::eigenDecompose_
protectedinherited
Vdouble bpp::AbstractSubstitutionModel::freq_
protectedinherited

The vector $\pi_e$ of equilibrium frequencies.

Definition at line 121 of file AbstractSubstitutionModel.h.

Referenced by bpp::AbstractSubstitutionModel::AbstractSubstitutionModel(), bpp::WordSubstitutionModel::completeMatrices(), bpp::Coala::computeEquilibriumFrequencies(), bpp::RE08::d2Pij_dt2(), bpp::RE08::dPij_dt(), bpp::DSO78::DSO78(), bpp::LG08::fireParameterChanged(), bpp::DSO78::fireParameterChanged(), bpp::WAG01::fireParameterChanged(), bpp::JTT92::fireParameterChanged(), bpp::UserProteinSubstitutionModel::fireParameterChanged(), bpp::JCprot::fireParameterChanged(), bpp::RE08::getd2Pij_dt2(), bpp::RE08::getdPij_dt(), bpp::AbstractSubstitutionModel::getFrequencies(), bpp::AbstractSubstitutionModel::getFrequencies_(), bpp::RN95s::getPij_t(), bpp::RN95::getPij_t(), bpp::RE08::getPij_t(), bpp::AbstractSubstitutionModel::getScale(), bpp::JCprot::JCprot(), bpp::JTT92::JTT92(), bpp::LG08::LG08(), bpp::AbstractSubstitutionModel::operator=(), bpp::RN95s::Pij_t(), bpp::RN95::Pij_t(), bpp::RE08::Pij_t(), bpp::RE08::RE08(), bpp::UserProteinSubstitutionModel::readFromFile(), bpp::AbstractSubstitutionModel::setFreq(), bpp::LG08::setFreqFromData(), bpp::DSO78::setFreqFromData(), bpp::JTT92::setFreqFromData(), bpp::WAG01::setFreqFromData(), bpp::UserProteinSubstitutionModel::setFreqFromData(), bpp::JCprot::setFreqFromData(), bpp::L95::updateMatrices(), bpp::gBGC::updateMatrices(), bpp::SSR::updateMatrices(), bpp::MixtureOfASubstitutionModel::updateMatrices(), bpp::RN95s::updateMatrices(), bpp::BinarySubstitutionModel::updateMatrices(), bpp::YpR::updateMatrices(), bpp::GTR::updateMatrices(), bpp::AbstractWordSubstitutionModel::updateMatrices(), bpp::TN93::updateMatrices(), bpp::RN95::updateMatrices(), bpp::JCnuc::updateMatrices(), updateMatrices(), bpp::K80::updateMatrices(), bpp::T92::updateMatrices(), bpp::RE08::updateMatrices(), bpp::HKY85::updateMatrices(), bpp::F84::updateMatrices(), bpp::JCprot::updateMatrices(), bpp::AbstractSubstitutionModel::updateMatrices(), bpp::AbstractReversibleSubstitutionModel::updateMatrices(), bpp::UserProteinSubstitutionModel::UserProteinSubstitutionModel(), and bpp::WAG01::WAG01().

bool bpp::AbstractSubstitutionModel::isScalable_
protectedinherited
RowMatrix<double> bpp::AbstractSubstitutionModel::pijt_
mutableprotectedinherited
double bpp::AbstractSubstitutionModel::rate_
protectedinherited

The rate of the model (default: 1). The generator (and all its vectorial components) is independent of the rate, since it should be normalized.

Definition at line 111 of file AbstractSubstitutionModel.h.

Referenced by bpp::AbstractSubstitutionModel::addRateParameter(), bpp::BinarySubstitutionModel::d2Pij_dt2(), bpp::RN95s::d2Pij_dt2(), bpp::JCnuc::d2Pij_dt2(), bpp::TN93::d2Pij_dt2(), bpp::RN95::d2Pij_dt2(), bpp::K80::d2Pij_dt2(), bpp::T92::d2Pij_dt2(), bpp::JCprot::d2Pij_dt2(), bpp::F84::d2Pij_dt2(), bpp::HKY85::d2Pij_dt2(), bpp::BinarySubstitutionModel::dPij_dt(), bpp::RN95s::dPij_dt(), bpp::JCnuc::dPij_dt(), bpp::TN93::dPij_dt(), bpp::RN95::dPij_dt(), bpp::K80::dPij_dt(), bpp::T92::dPij_dt(), bpp::JCprot::dPij_dt(), bpp::F84::dPij_dt(), bpp::HKY85::dPij_dt(), bpp::WordSubstitutionModel::getd2Pij_dt2(), bpp::BinarySubstitutionModel::getd2Pij_dt2(), bpp::RN95s::getd2Pij_dt2(), bpp::JCnuc::getd2Pij_dt2(), bpp::TN93::getd2Pij_dt2(), bpp::RN95::getd2Pij_dt2(), bpp::K80::getd2Pij_dt2(), bpp::T92::getd2Pij_dt2(), bpp::JCprot::getd2Pij_dt2(), bpp::F84::getd2Pij_dt2(), bpp::HKY85::getd2Pij_dt2(), bpp::AbstractSubstitutionModel::getd2Pij_dt2(), bpp::WordSubstitutionModel::getdPij_dt(), bpp::BinarySubstitutionModel::getdPij_dt(), bpp::RN95s::getdPij_dt(), bpp::JCnuc::getdPij_dt(), bpp::TN93::getdPij_dt(), bpp::RN95::getdPij_dt(), bpp::K80::getdPij_dt(), bpp::T92::getdPij_dt(), bpp::JCprot::getdPij_dt(), bpp::F84::getdPij_dt(), bpp::HKY85::getdPij_dt(), bpp::AbstractSubstitutionModel::getdPij_dt(), bpp::WordSubstitutionModel::getPij_t(), bpp::BinarySubstitutionModel::getPij_t(), bpp::RN95s::getPij_t(), bpp::JCnuc::getPij_t(), bpp::TN93::getPij_t(), bpp::RN95::getPij_t(), bpp::K80::getPij_t(), bpp::T92::getPij_t(), bpp::JCprot::getPij_t(), bpp::F84::getPij_t(), bpp::HKY85::getPij_t(), bpp::AbstractSubstitutionModel::getPij_t(), bpp::AbstractSubstitutionModel::getRate(), bpp::AbstractMixedSubstitutionModel::normalizeVRates(), bpp::AbstractSubstitutionModel::operator=(), bpp::BinarySubstitutionModel::Pij_t(), bpp::RN95s::Pij_t(), bpp::JCnuc::Pij_t(), bpp::TN93::Pij_t(), bpp::RN95::Pij_t(), bpp::K80::Pij_t(), bpp::T92::Pij_t(), bpp::JCprot::Pij_t(), bpp::F84::Pij_t(), bpp::HKY85::Pij_t(), bpp::AbstractMixedSubstitutionModel::setRate(), bpp::AbstractSubstitutionModel::setRate(), bpp::BinarySubstitutionModel::updateMatrices(), and updateMatrices().

std::shared_ptr<const StateMap> bpp::AbstractSubstitutionModel::stateMap_
protectedinherited
std::vector<double> bpp::AbstractMixedSubstitutionModel::vRates_
protectedinherited

vector of the rates of the models.

For the computation of the transition probabilities, the rates are included in the submodels while updating the mixture, so there is no need to multiply here the transition times with the rates.

The mean (on the distribution of the models) of the elements of this vector equals the overall rate of the mixture model, that is rate_;

Definition at line 97 of file AbstractMixedSubstitutionModel.h.

Referenced by bpp::AbstractMixedSubstitutionModel::AbstractMixedSubstitutionModel(), bpp::AbstractMixedSubstitutionModel::getVRates(), bpp::MixtureOfASubstitutionModel::MixtureOfASubstitutionModel(), MixtureOfSubstitutionModels(), bpp::AbstractMixedSubstitutionModel::normalizeVRates(), bpp::AbstractMixedSubstitutionModel::operator=(), bpp::AbstractMixedSubstitutionModel::setRate(), bpp::AbstractMixedSubstitutionModel::setVRates(), setVRates(), and updateMatrices().


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