bpp-phyl3  3.0.0
bpp::AbstractSubstitutionModel Class Referenceabstract

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

+ Inheritance diagram for bpp::AbstractSubstitutionModel:
+ Collaboration diagram for bpp::AbstractSubstitutionModel:

Public Member Functions

 AbstractSubstitutionModel (std::shared_ptr< const Alphabet > alpha, std::shared_ptr< const StateMapInterface > stateMap, const std::string &prefix)
 
 AbstractSubstitutionModel (const AbstractSubstitutionModel &model)
 
AbstractSubstitutionModeloperator= (const AbstractSubstitutionModel &model)
 
virtual ~AbstractSubstitutionModel ()
 
bool computeFrequencies () const
 
void computeFrequencies (bool yn)
 
const Matrix< double > & generator () const
 
const Matrix< double > & exchangeabilityMatrix () const
 
const Matrix< double > & getPij_t (double t) const
 
const Matrix< double > & getdPij_dt (double t) const
 
const Matrix< double > & getd2Pij_dt2 (double t) const
 
double Sij (size_t i, size_t j) const
 
const VdoublegetEigenValues () const
 
const VdoublegetIEigenValues () const
 
bool isDiagonalizable () const
 
bool isNonSingular () const
 
const Matrix< double > & getRowLeftEigenVectors () const
 
const Matrix< double > & getColumnRightEigenVectors () const
 
virtual double Qij (size_t i, size_t j) const
 A method for computing all necessary matrices. 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...
 
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
 returns if model is scalable More...
 
double getScale () const
 return scale More...
 
void setScale (double scale)
 Multiplies the current generator by the given scale. More...
 
void normalize ()
 normalize the generator More...
 
void setDiagonal ()
 set the diagonal of the generator such that sum on each line equals 0. More...
 
const Alphabetalphabet () const override
 
std::shared_ptr< const AlphabetgetAlphabet () const override
 
const StateMapInterfacestateMap () const override
 
std::shared_ptr< const StateMapInterfacegetStateMap () const override
 
size_t getNumberOfStates () const override
 Get the number of states. More...
 
const std::vector< int > & getAlphabetStates () const override
 
std::string getAlphabetStateAsChar (size_t index) const override
 
int getAlphabetStateAsInt (size_t index) const override
 
std::vector< size_t > getModelStates (int code) const override
 Get the state in the model corresponding to a particular state in the alphabet. More...
 
std::vector< size_t > getModelStates (const std::string &code) const override
 Get the state in the model corresponding to a particular state in the alphabet. More...
 
const VdoublegetFrequencies () const override
 
virtual double freq (size_t i) const override
 
virtual double Pij_t (size_t i, size_t j, double t) const override
 
virtual double dPij_dt (size_t i, size_t j, double t) const override
 
virtual double d2Pij_dt2 (size_t i, size_t j, double t) const override
 
double getInitValue (size_t i, int state) const override
 
void setFreqFromData (const SequenceDataInterface &data, double pseudoCount=0) override
 Set equilibrium frequencies equal to the frequencies estimated from the data. More...
 
virtual void setFreq (std::map< int, double > &freqs) override
 Set equilibrium frequencies. More...
 
const FrequencySetInterfacefrequencySet () const override
 
virtual void fireParameterChanged (const ParameterList &parameters) override
 Tells the model that a parameter value has changed. More...
 
void addRateParameter () override
 add a "rate" parameter to the model, that handles the overall rate of the process. More...
 
void setVerboseLevel (short level)
 
short verboseLevel () const
 
virtual double getRate () const override
 The rate of the substitution process. More...
 
virtual void setRate (double rate) override
 Set the rate of the model (must be positive). More...
 
const Eigen::VectorXd & Lik_t (const Eigen::VectorXd &values, double t) const override
 
const Eigen::VectorXd & dLik_dt (const Eigen::VectorXd &values, double t) const override
 
const Eigen::VectorXd & d2Lik_dt2 (const Eigen::VectorXd &values, double t) const override
 
TransitionModelInterfaceclone () const =0
 
virtual Clonableclone () const=0
 
virtual std::string getName () const =0
 Get the name of the model. 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 override
 
virtual const ParameterListgetParameters () const=0
 
const ParameterListgetParameters () const override
 
virtual const Parameterparameter (const std::string &name) const=0
 
const Parameterparameter (const std::string &name) const override
 
virtual double getParameterValue (const std::string &name) const=0
 
double getParameterValue (const std::string &name) const override
 
virtual void setAllParametersValues (const ParameterList &parameters)=0
 
void setAllParametersValues (const ParameterList &parameters) override
 
virtual void setParameterValue (const std::string &name, double value)=0
 
void setParameterValue (const std::string &name, double value) override
 
virtual void setParametersValues (const ParameterList &parameters)=0
 
void setParametersValues (const ParameterList &parameters) override
 
virtual bool matchParametersValues (const ParameterList &parameters)=0
 
bool matchParametersValues (const ParameterList &parameters) override
 
virtual void removeConstraint (const std::string &name)=0
 
void removeConstraint (const std::string &name) override
 
virtual void setConstraint (const std::string &name, std::shared_ptr< ConstraintInterface > constraint)=0
 
void setConstraint (const std::string &name, std::shared_ptr< ConstraintInterface > constraint) override
 
virtual size_t getNumberOfParameters () const=0
 
size_t getNumberOfParameters () const override
 
virtual void setNamespace (const std::string &prefix)=0
 
void setNamespace (const std::string &prefix)
 
virtual std::string getNamespace () const=0
 
std::string getNamespace () const override
 
virtual std::string getParameterNameWithoutNamespace (const std::string &name) const=0
 
std::string getParameterNameWithoutNamespace (const std::string &name) const override
 
bool hasIndependentParameter (const std::string &name) const
 
ParameterList getAliasedParameters (const ParameterList &pl) const
 
ParameterList getFromParameters (const ParameterList &pl) const
 
std::string getFrom (const std::string &name) const
 
const std::shared_ptr< Parameter > & getParameter (const std::string &name) const
 
SubstitutionModelInterfaceclone () const override=0
 

Protected Member Functions

virtual void updateMatrices_ ()
 Diagonalize the $Q$ matrix, and fill the eigenValues_, iEigenValues_, leftEigenVectors_ and rightEigenVectors_ matrices. More...
 
VdoublegetFrequencies_ () override
 
virtual ParameterListgetParameters_ ()=0
 
ParameterListgetParameters_ () override
 
const std::shared_ptr< Parameter > & getParameter (size_t i) const
 
std::shared_ptr< Parameter > & getParameter (size_t i)
 
void addParameter_ (Parameter *parameter)
 
void addParameters_ (const ParameterList &parameters)
 
void shareParameter_ (const std::shared_ptr< Parameter > &parameter)
 
void shareParameters_ (const ParameterList &parameters)
 
void includeParameters_ (const ParameterList &parameters)
 
void deleteParameter_ (size_t index)
 
void deleteParameter_ (std::string &name)
 
void deleteParameters_ (const std::vector< std::string > &names)
 
void resetParameters_ ()
 
ParametergetParameter_ (const std::string &name)
 
ParametergetParameter_ (size_t index)
 
const ParametergetParameter_ (size_t index) const
 
ParametergetParameterWithNamespace_ (const std::string &name)
 
const ParametergetParameterWithNamespace_ (const std::string &name) const
 

Protected Attributes

bool isScalable_
 If the model is scalable (ie generator can be normalized automatically). More...
 
RowMatrix< double > generator_
 The generator matrix $Q$ of the model. 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 symmetric. More...
 
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...
 
std::shared_ptr< const Alphabetalphabet_
 The alphabet relevant to this model. More...
 
std::shared_ptr< const StateMapInterfacestateMap_
 The map of model states with alphabet states. More...
 
size_t size_
 The number of states. 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...
 
Vdouble freq_
 The vector $\pi_e$ of equilibrium frequencies. More...
 
RowMatrix< double > pijt_
 These ones are for bookkeeping: More...
 
RowMatrix< double > dpijt_
 
RowMatrix< double > d2pijt_
 
short verboseLevel_
 

Private Attributes

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

Friends

class OneChangeRegisterTransitionModel
 

Detailed Description

Definition at line 284 of file AbstractSubstitutionModel.h.

Constructor & Destructor Documentation

◆ AbstractSubstitutionModel() [1/2]

AbstractSubstitutionModel::AbstractSubstitutionModel ( std::shared_ptr< const Alphabet alpha,
std::shared_ptr< const StateMapInterface stateMap,
const std::string &  prefix 
)

Definition at line 112 of file AbstractSubstitutionModel.cpp.

◆ AbstractSubstitutionModel() [2/2]

bpp::AbstractSubstitutionModel::AbstractSubstitutionModel ( const AbstractSubstitutionModel model)
inline

Definition at line 365 of file AbstractSubstitutionModel.h.

◆ ~AbstractSubstitutionModel()

virtual bpp::AbstractSubstitutionModel::~AbstractSubstitutionModel ( )
inlinevirtual

Definition at line 401 of file AbstractSubstitutionModel.h.

Member Function Documentation

◆ addRateParameter()

void AbstractTransitionModel::addRateParameter ( )
overridevirtualinherited

add a "rate" parameter to the model, that handles the overall rate of the process.

Implements bpp::BranchModelInterface.

Reimplemented in bpp::RegisterRatesSubstitutionModel.

Definition at line 63 of file AbstractSubstitutionModel.cpp.

References bpp::AbstractParameterAliasable::addParameter_(), bpp::AbstractParameterAliasable::getNamespace(), bpp::Parameter::R_PLUS_STAR, and bpp::AbstractTransitionModel::rate_.

◆ alphabet()

const Alphabet& bpp::AbstractTransitionModel::alphabet ( ) const
inlineoverridevirtualinherited

◆ clone() [1/2]

TransitionModelInterface* bpp::TransitionModelInterface::clone ( ) const
pure virtualinherited

Implements bpp::BranchModelInterface.

Implemented in bpp::ReversibleSubstitutionModelInterface, bpp::SubstitutionModelInterface, bpp::AbstractReversibleProteinSubstitutionModel, bpp::AbstractProteinSubstitutionModel, bpp::ProteinReversibleSubstitutionModelInterface, bpp::ProteinSubstitutionModelInterface, bpp::AbstractReversibleNucleotideSubstitutionModel, bpp::AbstractNucleotideSubstitutionModel, bpp::NucleotideReversibleSubstitutionModelInterface, bpp::NucleotideSubstitutionModelInterface, bpp::MixedTransitionModelInterface, bpp::MarkovModulatedSubstitutionModel, bpp::YNGP_M, bpp::CodonReversibleSubstitutionModelInterface, bpp::CodonSubstitutionModelInterface, bpp::AbstractKroneckerCodonSubstitutionModel, bpp::AbstractDFPSubstitutionModel, bpp::AbstractCodonSubstitutionModel, bpp::AbstractReversibleSubstitutionModel, bpp::AbstractMixedTransitionModel, bpp::WordSubstitutionModel, bpp::TS98, bpp::RegisterRatesSubstitutionModel, bpp::RE08Codon, bpp::RE08Protein, bpp::RE08Nucleotide, 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::MixtureOfTransitionModels, bpp::MixtureOfSubstitutionModels, bpp::MixtureOfATransitionModel, bpp::MixtureOfASubstitutionModel, 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::TwoParameterBinarySubstitutionModel, bpp::LG10_EX_EHO, bpp::LG10_EX_EHO::EmbeddedModel, bpp::HKY85, bpp::F84, bpp::CodonDistanceSubstitutionModel, and bpp::BinarySubstitutionModel.

Referenced by bpp::AbstractSinglePhyloSubstitutionMapping::addModel(), bpp::MixtureOfATransitionModel::MixtureOfATransitionModel(), and bpp::AbstractMixedTransitionModel::operator=().

◆ clone() [2/2]

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

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]

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

Implements bpp::TransitionModelInterface.

Definition at line 406 of file AbstractSubstitutionModel.h.

References computeFreq_.

◆ d2Lik_dt2()

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

◆ d2Pij_dt2()

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

Implements bpp::TransitionModelInterface.

Reimplemented in bpp::RE08, bpp::JCprot, bpp::TN93, bpp::T92, bpp::K80, bpp::JCnuc, bpp::F81, bpp::EquiprobableSubstitutionModel, bpp::TwoParameterBinarySubstitutionModel, bpp::HKY85, bpp::F84, and bpp::BinarySubstitutionModel.

Definition at line 199 of file AbstractSubstitutionModel.h.

References bpp::AbstractTransitionModel::getd2Pij_dt2().

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

◆ dLik_dt()

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

◆ dPij_dt()

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

Implements bpp::TransitionModelInterface.

Reimplemented in bpp::RE08, bpp::JCprot, bpp::TN93, bpp::T92, bpp::K80, bpp::JCnuc, bpp::F81, bpp::EquiprobableSubstitutionModel, bpp::TwoParameterBinarySubstitutionModel, bpp::HKY85, bpp::F84, and bpp::BinarySubstitutionModel.

Definition at line 198 of file AbstractSubstitutionModel.h.

References bpp::AbstractTransitionModel::getdPij_dt().

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

◆ enableEigenDecomposition() [1/2]

◆ enableEigenDecomposition() [2/2]

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

Set if eigenValues and Vectors must be computed.

Implements bpp::SubstitutionModelInterface.

Definition at line 432 of file AbstractSubstitutionModel.h.

References eigenDecompose_.

◆ exchangeabilityMatrix()

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

Implements bpp::SubstitutionModelInterface.

Definition at line 410 of file AbstractSubstitutionModel.h.

References exchangeability_.

◆ fireParameterChanged()

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

Tells the model that a parameter value has changed.

This updates the matrices consequently.

Reimplemented from bpp::AbstractParameterAliasable.

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

Definition at line 217 of file AbstractSubstitutionModel.h.

References bpp::AbstractParameterAliasable::fireParameterChanged(), bpp::AbstractParameterAliasable::getNamespace(), bpp::ParameterList::getParameterValue(), bpp::ParameterList::hasParameter(), bpp::AbstractTransitionModel::rate_, bpp::ParameterList::size(), and bpp::AbstractTransitionModel::updateMatrices_().

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

◆ freq()

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

Implements bpp::TransitionModelInterface.

Definition at line 195 of file AbstractSubstitutionModel.h.

References bpp::AbstractTransitionModel::freq_.

Referenced by bpp::CodonSameAARateSubstitutionModel::compute_().

◆ frequencySet()

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

Implements bpp::BranchModelInterface.

Reimplemented in bpp::RegisterRatesSubstitutionModel, bpp::WAG01, bpp::UserProteinSubstitutionModel, bpp::LG08, bpp::JTT92, bpp::JCprot, bpp::DSO78, bpp::EquiprobableSubstitutionModel, bpp::D1WalkSubstitutionModel, and bpp::CodonDistanceFrequenciesSubstitutionModel.

Definition at line 207 of file AbstractSubstitutionModel.h.

◆ generator()

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

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

See also
Qij()

Implements bpp::SubstitutionModelInterface.

Definition at line 408 of file AbstractSubstitutionModel.h.

References generator_.

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

◆ getAlphabet()

std::shared_ptr<const Alphabet> bpp::AbstractTransitionModel::getAlphabet ( ) const
inlineoverridevirtualinherited

◆ getAlphabetStateAsChar()

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

Implements bpp::BranchModelInterface.

Reimplemented in bpp::RegisterRatesSubstitutionModel.

Definition at line 179 of file AbstractSubstitutionModel.h.

References bpp::AbstractTransitionModel::stateMap_.

◆ getAlphabetStateAsInt()

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

Implements bpp::BranchModelInterface.

Reimplemented in bpp::RegisterRatesSubstitutionModel.

Definition at line 181 of file AbstractSubstitutionModel.h.

References bpp::AbstractTransitionModel::stateMap_.

Referenced by bpp::AbstractCodonSubstitutionModel::completeMatrices_(), and bpp::AbstractTransitionModel::getInitValue().

◆ getAlphabetStates()

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

Implements bpp::BranchModelInterface.

Reimplemented in bpp::RegisterRatesSubstitutionModel.

Definition at line 177 of file AbstractSubstitutionModel.h.

References bpp::AbstractTransitionModel::stateMap_.

Referenced by bpp::DFP07::DFP07(), bpp::YNGP_M3::YNGP_M3(), and bpp::YNGP_M8::YNGP_M8().

◆ getColumnRightEigenVectors()

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

Implements bpp::SubstitutionModelInterface.

Definition at line 428 of file AbstractSubstitutionModel.h.

References rightEigenVectors_.

◆ getd2Pij_dt2()

◆ getdPij_dt()

◆ getEigenValues()

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

Implements bpp::SubstitutionModelInterface.

Definition at line 418 of file AbstractSubstitutionModel.h.

References eigenValues_.

◆ getFrequencies()

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

Implements bpp::TransitionModelInterface.

Definition at line 187 of file AbstractSubstitutionModel.h.

References bpp::AbstractTransitionModel::freq_.

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

◆ getFrequencies_()

◆ getIEigenValues()

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

Implements bpp::SubstitutionModelInterface.

Definition at line 420 of file AbstractSubstitutionModel.h.

References iEigenValues_.

◆ getInitValue()

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

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

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

Implements bpp::BranchModelInterface.

Reimplemented in bpp::RE08.

Definition at line 70 of file AbstractSubstitutionModel.cpp.

References bpp::AbstractTransitionModel::alphabet_, bpp::AbstractTransitionModel::getAlphabetStateAsInt(), and bpp::AbstractTransitionModel::size_.

◆ getModelStates() [1/2]

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

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

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

Implements bpp::BranchModelInterface.

Reimplemented in bpp::RegisterRatesSubstitutionModel.

Definition at line 185 of file AbstractSubstitutionModel.h.

References bpp::AbstractTransitionModel::stateMap_.

◆ getModelStates() [2/2]

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

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

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

Implements bpp::BranchModelInterface.

Reimplemented in bpp::RegisterRatesSubstitutionModel.

Definition at line 183 of file AbstractSubstitutionModel.h.

References bpp::AbstractTransitionModel::stateMap_.

◆ 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(), getPij_t(), bpp::BppOSubstitutionModelFormat::initialize_(), bpp::MixtureOfSubstitutionModels::MixtureOfSubstitutionModels(), bpp::MixtureOfATransitionModel::model(), bpp::LegacyPhylogeneticsApplicationTools::setSubstitutionModelParametersInitialValuesWithAliases(), bpp::ModelPath::toString(), and bpp::BppOSubstitutionModelFormat::write().

◆ getNumberOfStates()

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

◆ getPij_t()

◆ getRate()

double AbstractTransitionModel::getRate ( ) const
overridevirtualinherited

The rate of the substitution process.

Implements bpp::BranchModelInterface.

Reimplemented in bpp::RegisterRatesSubstitutionModel.

Definition at line 45 of file AbstractSubstitutionModel.cpp.

References bpp::AbstractTransitionModel::rate_.

◆ getRowLeftEigenVectors()

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

Implements bpp::SubstitutionModelInterface.

Definition at line 426 of file AbstractSubstitutionModel.h.

References leftEigenVectors_.

◆ getScale()

double AbstractSubstitutionModel::getScale ( ) const
virtual

◆ getStateMap()

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

Implements bpp::BranchModelInterface.

Reimplemented in bpp::RegisterRatesSubstitutionModel.

Definition at line 173 of file AbstractSubstitutionModel.h.

References bpp::AbstractTransitionModel::stateMap_.

Referenced by bpp::D1WalkSubstitutionModel::D1WalkSubstitutionModel(), and bpp::EquiprobableSubstitutionModel::EquiprobableSubstitutionModel().

◆ isDiagonalizable()

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

Implements bpp::SubstitutionModelInterface.

Definition at line 422 of file AbstractSubstitutionModel.h.

References isDiagonalizable_.

◆ isNonSingular()

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

Implements bpp::SubstitutionModelInterface.

Definition at line 424 of file AbstractSubstitutionModel.h.

References isNonSingular_.

◆ isScalable()

virtual bool bpp::AbstractSubstitutionModel::isScalable ( ) const
inlinevirtual

◆ Lik_t()

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

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

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

Implements bpp::TransitionModelInterface.

Definition at line 30 of file AbstractSubstitutionModel.h.

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

◆ normalize()

void AbstractSubstitutionModel::normalize ( )
virtual

◆ operator=()

◆ Pij_t()

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

Implements bpp::TransitionModelInterface.

Reimplemented in bpp::RE08, bpp::JCprot, bpp::TN93, bpp::T92, bpp::K80, bpp::JCnuc, bpp::F81, bpp::EquiprobableSubstitutionModel, bpp::TwoParameterBinarySubstitutionModel, bpp::HKY85, bpp::F84, and bpp::BinarySubstitutionModel.

Definition at line 197 of file AbstractSubstitutionModel.h.

References bpp::AbstractTransitionModel::getPij_t().

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

◆ Qij()

virtual double bpp::AbstractSubstitutionModel::Qij ( size_t  i,
size_t  j 
) const
inlinevirtual

A method for computing all necessary matrices.

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

Implements bpp::SubstitutionModelInterface.

Definition at line 430 of file AbstractSubstitutionModel.h.

References generator_.

◆ setDiagonal()

◆ setFreq()

◆ setFreqFromData()

void AbstractTransitionModel::setFreqFromData ( const SequenceDataInterface data,
double  pseudoCount = 0 
)
overridevirtualinherited

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

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

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

Implements bpp::TransitionModelInterface.

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

Definition at line 88 of file AbstractSubstitutionModel.cpp.

References bpp::SequenceContainerTools::getFrequencies(), and bpp::AbstractTransitionModel::setFreq().

◆ setRate()

void AbstractTransitionModel::setRate ( double  rate)
overridevirtualinherited

◆ setScalable()

void bpp::AbstractSubstitutionModel::setScalable ( bool  scalable)
inlinevirtual

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

Implements bpp::SubstitutionModelInterface.

Definition at line 461 of file AbstractSubstitutionModel.h.

References isScalable_.

◆ setScale()

◆ setVerboseLevel()

void bpp::AbstractTransitionModel::setVerboseLevel ( short  level)
inlineinherited

◆ Sij()

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

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

Implements bpp::SubstitutionModelInterface.

Definition at line 416 of file AbstractSubstitutionModel.h.

References exchangeability_.

◆ stateMap()

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

Implements bpp::BranchModelInterface.

Reimplemented in bpp::RegisterRatesSubstitutionModel.

Definition at line 171 of file AbstractSubstitutionModel.h.

References bpp::AbstractTransitionModel::stateMap_.

◆ updateMatrices_()

void AbstractSubstitutionModel::updateMatrices_ ( )
protectedvirtual

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

Implements bpp::AbstractTransitionModel.

Reimplemented in bpp::WordSubstitutionModel, bpp::RegisterRatesSubstitutionModel, bpp::RE08, bpp::JCprot, bpp::Coala, bpp::POMO, bpp::YpR_Gen, bpp::YpR_Sym, bpp::YpR, bpp::TN93, bpp::T92, bpp::SSR, bpp::RN95s, bpp::RN95, bpp::L95, bpp::K80, bpp::JCnuc, bpp::GTR, bpp::gBGC, bpp::F81, bpp::EquiprobableSubstitutionModel, bpp::D1WalkSubstitutionModel, bpp::AbstractDFPSubstitutionModel, bpp::AbstractCodonSubstitutionModel, bpp::AbstractReversibleSubstitutionModel, bpp::TwoParameterBinarySubstitutionModel, bpp::HKY85, bpp::F84, bpp::BinarySubstitutionModel, and bpp::AbstractWordSubstitutionModel.

Definition at line 136 of file AbstractSubstitutionModel.cpp.

References bpp::abs(), bpp::MatrixTools::add(), computeFrequencies(), bpp::ApplicationTools::displayMessage(), eigenValues_, enableEigenDecomposition(), bpp::AbstractTransitionModel::freq_, generator_, bpp::MatrixTools::getId(), bpp::EigenValue< class >::getImagEigenValues(), bpp::AbstractTransitionModel::getNumberOfStates(), bpp::EigenValue< class >::getRealEigenValues(), bpp::EigenValue< class >::getV(), iEigenValues_, bpp::MatrixTools::inv(), isDiagonalizable_, isNonSingular_, leftEigenVectors_, bpp::NumConstants::NANO(), normalize(), bpp::MatrixTools::pow(), RowMatrix< double >::resize(), rightEigenVectors_, setScale(), bpp::NumConstants::SMALL(), bpp::VectorTools::sum(), bpp::MatrixTools::Taylor(), bpp::NumConstants::TINY(), tmpMat_, bpp::AbstractTransitionModel::verboseLevel_, and vPowGen_.

Referenced by bpp::CodonSameAARateSubstitutionModel::CodonSameAARateSubstitutionModel(), bpp::CodonSameAARateSubstitutionModel::fireParameterChanged(), bpp::AbstractWordSubstitutionModel::updateMatrices_(), bpp::AbstractReversibleSubstitutionModel::updateMatrices_(), bpp::AbstractDFPSubstitutionModel::updateMatrices_(), bpp::L95::updateMatrices_(), bpp::POMO::updateMatrices_(), bpp::RE08::updateMatrices_(), and bpp::RegisterRatesSubstitutionModel::updateMatrices_().

◆ verboseLevel()

short bpp::AbstractTransitionModel::verboseLevel ( ) const
inlineinherited

Friends And Related Function Documentation

◆ OneChangeRegisterTransitionModel

friend class OneChangeRegisterTransitionModel
friend

Definition at line 497 of file AbstractSubstitutionModel.h.

Member Data Documentation

◆ alphabet_

◆ computeFreq_

bool bpp::AbstractSubstitutionModel::computeFreq_
protected

if the Frequencies must be computed from the generator

Definition at line 303 of file AbstractSubstitutionModel.h.

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

◆ d2pijt_

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

◆ dpijt_

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

◆ eigenDecompose_

bool bpp::AbstractSubstitutionModel::eigenDecompose_
protected

Tell if the eigen decomposition should be performed.

Definition at line 315 of file AbstractSubstitutionModel.h.

Referenced by enableEigenDecomposition(), and operator=().

◆ eigenValues_

◆ exchangeability_

◆ freq_

Vdouble bpp::AbstractTransitionModel::freq_
protectedinherited

The vector $\pi_e$ of equilibrium frequencies.

Definition at line 143 of file AbstractSubstitutionModel.h.

Referenced by bpp::AbstractReversibleSubstitutionModel::AbstractReversibleSubstitutionModel(), bpp::AbstractTransitionModel::AbstractTransitionModel(), bpp::WordSubstitutionModel::completeMatrices_(), bpp::Coala::computeEquilibriumFrequencies(), bpp::RE08::d2Pij_dt2(), bpp::RE08::dPij_dt(), bpp::DSO78::DSO78(), bpp::EquiprobableSubstitutionModel::EquiprobableSubstitutionModel(), bpp::EquiprobableSubstitutionModel::fireParameterChanged(), bpp::DSO78::fireParameterChanged(), bpp::JCprot::fireParameterChanged(), bpp::JTT92::fireParameterChanged(), bpp::LG08::fireParameterChanged(), bpp::UserProteinSubstitutionModel::fireParameterChanged(), bpp::WAG01::fireParameterChanged(), bpp::AbstractTransitionModel::freq(), bpp::RE08::getd2Pij_dt2(), bpp::RE08::getdPij_dt(), bpp::AbstractTransitionModel::getFrequencies(), bpp::AbstractTransitionModel::getFrequencies_(), bpp::RE08::getPij_t(), getScale(), bpp::JCprot::JCprot(), bpp::JTT92::JTT92(), bpp::LG08::LG08(), bpp::RE08::Pij_t(), bpp::RE08::RE08(), bpp::UserProteinSubstitutionModel::readFromFile(), bpp::D1WalkSubstitutionModel::setFreq(), bpp::EquiprobableSubstitutionModel::setFreq(), bpp::AbstractTransitionModel::setFreq(), bpp::DSO78::setFreqFromData(), bpp::JCprot::setFreqFromData(), bpp::JTT92::setFreqFromData(), bpp::LG08::setFreqFromData(), bpp::UserProteinSubstitutionModel::setFreqFromData(), bpp::WAG01::setFreqFromData(), updateMatrices_(), bpp::AbstractWordSubstitutionModel::updateMatrices_(), bpp::BinarySubstitutionModel::updateMatrices_(), bpp::F84::updateMatrices_(), bpp::HKY85::updateMatrices_(), bpp::TwoParameterBinarySubstitutionModel::updateMatrices_(), bpp::AbstractReversibleSubstitutionModel::updateMatrices_(), bpp::D1WalkSubstitutionModel::updateMatrices_(), bpp::EquiprobableSubstitutionModel::updateMatrices_(), bpp::MixtureOfATransitionModel::updateMatrices_(), bpp::MixtureOfTransitionModels::updateMatrices_(), bpp::F81::updateMatrices_(), bpp::gBGC::updateMatrices_(), bpp::GTR::updateMatrices_(), bpp::JCnuc::updateMatrices_(), bpp::K80::updateMatrices_(), bpp::L95::updateMatrices_(), bpp::RN95::updateMatrices_(), bpp::RN95s::updateMatrices_(), bpp::SSR::updateMatrices_(), bpp::T92::updateMatrices_(), bpp::TN93::updateMatrices_(), bpp::POMO::updateMatrices_(), bpp::JCprot::updateMatrices_(), bpp::RE08::updateMatrices_(), bpp::YpR::updateMatrices_(), bpp::UserProteinSubstitutionModel::UserProteinSubstitutionModel(), and bpp::WAG01::WAG01().

◆ generator_

◆ iEigenValues_

Vdouble bpp::AbstractSubstitutionModel::iEigenValues_
protected

The vector of the imaginary part of the eigen values.

Definition at line 325 of file AbstractSubstitutionModel.h.

Referenced by getd2Pij_dt2(), getdPij_dt(), getIEigenValues(), getPij_t(), operator=(), setScale(), updateMatrices_(), bpp::gBGC::updateMatrices_(), and bpp::YpR::updateMatrices_().

◆ isDiagonalizable_

◆ isNonSingular_

◆ isScalable_

bool bpp::AbstractSubstitutionModel::isScalable_
protected

If the model is scalable (ie generator can be normalized automatically).

Definition at line 293 of file AbstractSubstitutionModel.h.

Referenced by isScalable(), normalize(), operator=(), bpp::RegisterRatesSubstitutionModel::RegisterRatesSubstitutionModel(), setScalable(), and setScale().

◆ leftEigenVectors_

◆ lik_

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

◆ pijt_

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

These ones are for bookkeeping:

Definition at line 148 of file AbstractSubstitutionModel.h.

Referenced by bpp::WordSubstitutionModel::getPij_t(), getPij_t(), and bpp::AbstractMixedTransitionModel::getPij_t().

◆ rate_

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

Referenced by bpp::AbstractTransitionModel::addRateParameter(), bpp::BinarySubstitutionModel::d2Pij_dt2(), bpp::F84::d2Pij_dt2(), bpp::HKY85::d2Pij_dt2(), bpp::TwoParameterBinarySubstitutionModel::d2Pij_dt2(), bpp::EquiprobableSubstitutionModel::d2Pij_dt2(), bpp::F81::d2Pij_dt2(), bpp::JCnuc::d2Pij_dt2(), bpp::K80::d2Pij_dt2(), bpp::T92::d2Pij_dt2(), bpp::TN93::d2Pij_dt2(), bpp::JCprot::d2Pij_dt2(), bpp::BinarySubstitutionModel::dPij_dt(), bpp::F84::dPij_dt(), bpp::HKY85::dPij_dt(), bpp::TwoParameterBinarySubstitutionModel::dPij_dt(), bpp::EquiprobableSubstitutionModel::dPij_dt(), bpp::F81::dPij_dt(), bpp::JCnuc::dPij_dt(), bpp::K80::dPij_dt(), bpp::T92::dPij_dt(), bpp::TN93::dPij_dt(), bpp::JCprot::dPij_dt(), bpp::AbstractTransitionModel::fireParameterChanged(), bpp::BinarySubstitutionModel::getd2Pij_dt2(), bpp::F84::getd2Pij_dt2(), bpp::HKY85::getd2Pij_dt2(), bpp::TwoParameterBinarySubstitutionModel::getd2Pij_dt2(), bpp::EquiprobableSubstitutionModel::getd2Pij_dt2(), bpp::F81::getd2Pij_dt2(), bpp::JCnuc::getd2Pij_dt2(), bpp::K80::getd2Pij_dt2(), bpp::T92::getd2Pij_dt2(), bpp::TN93::getd2Pij_dt2(), bpp::JCprot::getd2Pij_dt2(), bpp::WordSubstitutionModel::getd2Pij_dt2(), getd2Pij_dt2(), bpp::BinarySubstitutionModel::getdPij_dt(), bpp::F84::getdPij_dt(), bpp::HKY85::getdPij_dt(), bpp::TwoParameterBinarySubstitutionModel::getdPij_dt(), bpp::EquiprobableSubstitutionModel::getdPij_dt(), bpp::F81::getdPij_dt(), bpp::JCnuc::getdPij_dt(), bpp::K80::getdPij_dt(), bpp::T92::getdPij_dt(), bpp::TN93::getdPij_dt(), bpp::JCprot::getdPij_dt(), bpp::WordSubstitutionModel::getdPij_dt(), getdPij_dt(), bpp::BinarySubstitutionModel::getPij_t(), bpp::F84::getPij_t(), bpp::HKY85::getPij_t(), bpp::TwoParameterBinarySubstitutionModel::getPij_t(), bpp::EquiprobableSubstitutionModel::getPij_t(), bpp::F81::getPij_t(), bpp::JCnuc::getPij_t(), bpp::K80::getPij_t(), bpp::T92::getPij_t(), bpp::TN93::getPij_t(), bpp::JCprot::getPij_t(), bpp::WordSubstitutionModel::getPij_t(), getPij_t(), bpp::AbstractTransitionModel::getRate(), bpp::AbstractMixedTransitionModel::normalizeVRates(), bpp::BinarySubstitutionModel::Pij_t(), bpp::F84::Pij_t(), bpp::HKY85::Pij_t(), bpp::TwoParameterBinarySubstitutionModel::Pij_t(), bpp::EquiprobableSubstitutionModel::Pij_t(), bpp::F81::Pij_t(), bpp::JCnuc::Pij_t(), bpp::K80::Pij_t(), bpp::T92::Pij_t(), bpp::TN93::Pij_t(), bpp::JCprot::Pij_t(), bpp::AbstractMixedTransitionModel::setRate(), bpp::AbstractTransitionModel::setRate(), bpp::BinarySubstitutionModel::updateMatrices_(), bpp::TwoParameterBinarySubstitutionModel::updateMatrices_(), and bpp::MixtureOfTransitionModels::updateMatrices_().

◆ rightEigenVectors_

◆ size_

◆ stateMap_

◆ tmpMat_

RowMatrix<double> bpp::AbstractSubstitutionModel::tmpMat_
mutableprotected

◆ verboseLevel_

short bpp::AbstractTransitionModel::verboseLevel_
protectedinherited

◆ vPowGen_

std::vector< RowMatrix<double> > bpp::AbstractSubstitutionModel::vPowGen_
protected

vector of the powers of generator_ for Taylor development (if rightEigenVectors_ is singular).

Definition at line 352 of file AbstractSubstitutionModel.h.

Referenced by getd2Pij_dt2(), getdPij_dt(), getPij_t(), operator=(), updateMatrices_(), bpp::gBGC::updateMatrices_(), bpp::RN95::updateMatrices_(), bpp::RN95s::updateMatrices_(), and bpp::YpR::updateMatrices_().


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