bpp-phyl3
3.0.0
|
Partial implementation of the Markov-modulated class of substitution models. More...
#include <Bpp/Phyl/Model/MarkovModulatedSubstitutionModel.h>
Public Member Functions | |
MarkovModulatedSubstitutionModel (std::unique_ptr< ReversibleSubstitutionModelInterface > model, unsigned int nbRates, bool normalizeRateChanges, const std::string &prefix) | |
Build a new MarkovModulatedSubstitutionModel object. More... | |
MarkovModulatedSubstitutionModel (const MarkovModulatedSubstitutionModel &model) | |
MarkovModulatedSubstitutionModel & | operator= (const MarkovModulatedSubstitutionModel &model) |
virtual | ~MarkovModulatedSubstitutionModel () |
MarkovModulatedSubstitutionModel * | clone () const override=0 |
const Alphabet & | alphabet () const override |
std::shared_ptr< const Alphabet > | getAlphabet () const override |
size_t | getNumberOfStates () const override |
Get the number of states. More... | |
const StateMapInterface & | stateMap () const override |
std::shared_ptr< const StateMapInterface > | getStateMap () const override |
const std::vector< int > & | getAlphabetStates () const override |
std::string | getAlphabetStateAsChar (size_t index) const override |
int | getAlphabetStateAsInt (size_t index) const override |
std::vector< size_t > | getModelStates (int code) const override |
Get the state in the model corresponding to a particular state in the alphabet. More... | |
std::vector< size_t > | getModelStates (const std::string &code) const override |
Get the state in the model corresponding to a particular state in the alphabet. More... | |
const Vdouble & | getFrequencies () const override |
const Matrix< double > & | exchangeabilityMatrix () const override |
const Matrix< double > & | generator () const override |
const Matrix< double > & | getPij_t (double t) const override |
const Matrix< double > & | getdPij_dt (double t) const override |
const Matrix< double > & | getd2Pij_dt2 (double t) const override |
const Vdouble & | getEigenValues () const override |
const Vdouble & | getIEigenValues () const override |
bool | isDiagonalizable () const override |
bool | isNonSingular () const override |
const Matrix< double > & | getRowLeftEigenVectors () const override |
const Matrix< double > & | getColumnRightEigenVectors () const override |
double | freq (size_t i) const override |
double | Sij (size_t i, size_t j) const override |
double | Qij (size_t i, size_t j) const override |
A method for computing all necessary matrices. More... | |
double | Pij_t (size_t i, size_t j, double t) const override |
double | dPij_dt (size_t i, size_t j, double t) const override |
double | d2Pij_dt2 (size_t i, size_t j, double t) const override |
double | getInitValue (size_t i, int state) const override |
void | setFreqFromData (const SequenceDataInterface &data, double pseudoCount=0) override |
Set equilibrium frequencies equal to the frequencies estimated from the data. More... | |
void | setFreq (std::map< int, double > &frequencies) override |
Set equilibrium frequencies. More... | |
const FrequencySetInterface & | frequencySet () const override |
const ReversibleSubstitutionModelInterface & | nestedModel () const |
size_t | getRate (size_t i) const |
Get the rate category corresponding to a particular state in the compound model. More... | |
double | getRate () const override |
Get the rate. More... | |
void | setRate (double rate) override |
Set the rate of the model (must be positive). More... | |
bool | isScalable () const override |
returns if model is scalable More... | |
void | setScalable (bool scalable) override |
sets if model is scalable, ie scale can be changed. Default : true, set to false to avoid normalization for example. More... | |
void | normalize () override |
Normalize the generator. More... | |
void | setDiagonal () override |
set the diagonal of the generator such that sum on each line equals 0. More... | |
double | getScale () const override |
Get the scalar product of diagonal elements of the generator and the frequencies vector. If the generator is normalized, then scale=1. Otherwise each element must be multiplied by 1/scale. More... | |
void | setScale (double scale) override |
Multiplies the current generator by the given scale. More... | |
void | enableEigenDecomposition (bool yn) override |
Set if eigenValues and Vectors must be computed. More... | |
bool | enableEigenDecomposition () override |
Tell if eigenValues and Vectors must be computed. More... | |
bool | computeFrequencies () const override |
void | computeFrequencies (bool yn) override |
virtual void | fireParameterChanged (const ParameterList ¶meters) override |
Tells the model that a parameter value has changed. More... | |
void | setNamespace (const std::string &prefix) override |
virtual const Eigen::VectorXd & | Lik_t (const Eigen::VectorXd &values, double t) const =0 |
virtual const Eigen::VectorXd & | dLik_dt (const Eigen::VectorXd &values, double t) const =0 |
virtual const Eigen::VectorXd & | d2Lik_dt2 (const Eigen::VectorXd &values, double t) const =0 |
virtual std::string | getName () const =0 |
Get the name of the model. More... | |
virtual void | addRateParameter ()=0 |
virtual size_t | getNumberOfIndependentParameters () const=0 |
virtual void | aliasParameters (const std::string &p1, const std::string &p2)=0 |
virtual void | aliasParameters (std::map< std::string, std::string > &unparsedParams, bool verbose)=0 |
virtual void | unaliasParameters (const std::string &p1, const std::string &p2)=0 |
virtual const ParameterList & | getIndependentParameters () const=0 |
virtual std::vector< std::string > | getAlias (const std::string &name) const=0 |
virtual std::map< std::string, std::string > | getAliases () const=0 |
virtual bool | hasParameter (const std::string &name) const=0 |
virtual const ParameterList & | getParameters () const=0 |
virtual const Parameter & | parameter (const std::string &name) const=0 |
virtual double | getParameterValue (const std::string &name) const=0 |
virtual void | setAllParametersValues (const ParameterList ¶meters)=0 |
virtual void | setParameterValue (const std::string &name, double value)=0 |
virtual void | setParametersValues (const ParameterList ¶meters)=0 |
virtual bool | matchParametersValues (const ParameterList ¶meters)=0 |
virtual void | removeConstraint (const std::string &name)=0 |
virtual void | setConstraint (const std::string &name, std::shared_ptr< ConstraintInterface > constraint)=0 |
virtual size_t | getNumberOfParameters () const=0 |
virtual std::string | getNamespace () const=0 |
virtual std::string | getParameterNameWithoutNamespace (const std::string &name) const=0 |
bool | hasIndependentParameter (const std::string &name) const |
const ParameterList & | getIndependentParameters () const |
size_t | getNumberOfIndependentParameters () const |
void | aliasParameters (const std::string &p1, const std::string &p2) |
void | aliasParameters (std::map< std::string, std::string > &unparsedParams, bool verbose) |
void | unaliasParameters (const std::string &p1, const std::string &p2) |
ParameterList | getAliasedParameters (const ParameterList &pl) const |
ParameterList | getFromParameters (const ParameterList &pl) const |
virtual std::vector< std::string > | getAlias (const std::string &name) const |
virtual std::map< std::string, std::string > | getAliases () const |
std::string | getFrom (const std::string &name) const |
bool | hasParameter (const std::string &name) const override |
const ParameterList & | getParameters () const override |
const Parameter & | parameter (const std::string &name) const override |
const std::shared_ptr< Parameter > & | getParameter (const std::string &name) const |
double | getParameterValue (const std::string &name) const override |
void | setAllParametersValues (const ParameterList ¶meters) override |
void | setParameterValue (const std::string &name, double value) override |
void | setParametersValues (const ParameterList ¶meters) override |
bool | matchParametersValues (const ParameterList ¶meters) override |
void | removeConstraint (const std::string &name) override |
void | setConstraint (const std::string &name, std::shared_ptr< ConstraintInterface > constraint) override |
size_t | getNumberOfParameters () const override |
std::string | getNamespace () const override |
std::string | getParameterNameWithoutNamespace (const std::string &name) const override |
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 |
Protected Member Functions | |
virtual void | updateMatrices_ () |
virtual void | updateRatesModel_ ()=0 |
Update the rates vector, generator and equilibrium frequencies. More... | |
Vdouble & | getFrequencies_ () override |
virtual ParameterList & | getParameters_ ()=0 |
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 ¶meters) |
void | shareParameter_ (const std::shared_ptr< Parameter > ¶meter) |
void | shareParameters_ (const ParameterList ¶meters) |
void | includeParameters_ (const ParameterList ¶meters) |
void | deleteParameter_ (size_t index) |
void | deleteParameter_ (std::string &name) |
void | deleteParameters_ (const std::vector< std::string > &names) |
void | resetParameters_ () |
Parameter & | getParameter_ (const std::string &name) |
Parameter & | getParameter_ (size_t index) |
const Parameter & | getParameter_ (size_t index) const |
Parameter & | getParameterWithNamespace_ (const std::string &name) |
const Parameter & | getParameterWithNamespace_ (const std::string &name) const |
ParameterList & | getParameters_ () override |
Protected Attributes | |
std::unique_ptr< ReversibleSubstitutionModelInterface > | model_ |
std::shared_ptr< const MarkovModulatedStateMap > | stateMap_ |
size_t | nbStates_ |
size_t | nbRates_ |
RowMatrix< double > | ratesGenerator_ |
RowMatrix< double > | generator_ |
The generator matrix of the model. More... | |
RowMatrix< double > | exchangeability_ |
The exchangeability matrix of the model. More... | |
RowMatrix< double > | leftEigenVectors_ |
The matrix made of left eigen vectors (by row). More... | |
RowMatrix< double > | rightEigenVectors_ |
The matrix made of right eigen vectors (by column). More... | |
Vdouble | eigenValues_ |
The vector of real parts of eigen values. More... | |
Vdouble | iEigenValues_ |
The vector of imaginary parts of the eigen values (zero in case of reversible pmodel). More... | |
bool | eigenDecompose_ |
Tell if the eigen decomposition should be performed. More... | |
bool | compFreq_ |
Tell if the equilibrium frequencies should be computed from the generator. More... | |
RowMatrix< double > | pijt_ |
These ones are for bookkeeping: More... | |
RowMatrix< double > | dpijt_ |
RowMatrix< double > | d2pijt_ |
Vdouble | freq_ |
The vector of equilibrium frequencies. More... | |
bool | normalizeRateChanges_ |
std::string | nestedPrefix_ |
Rate generator. | |
These variables must be initialized in the constructor of the derived class. | |
RowMatrix< double > | rates_ |
RowMatrix< double > | ratesExchangeability_ |
Vdouble | ratesFreq_ |
Private Attributes | |
ParameterList | independentParameters_ |
std::map< std::string, std::shared_ptr< AliasParameterListener > > | aliasListenersRegister_ |
ParameterList | parameters_ |
std::string | prefix_ |
Eigen::VectorXd | lik_ |
Partial implementation of the Markov-modulated class of substitution models.
This class wraps a substitution model and provide a matrix describing rate changes. The rate matrix must be initialized by derived classes of this class. Using these matrices, the diagonalization procedure of Galtier and Jean-Marie is used.
Such models can be described using two matrices: a substitution matrix, , with size , which is a "standard" substitution model of any alphabet type, and a rate matrix of size . The generator of the markov-modulated model, can be written using Kronecker matrix operands.
where is the diagonal matrix of all rates, and is the identity matrix of size .
This generator is normalized so that branch lengths are measured in unit of mean number of substitutions per site, where susbstitution here means "change of alphabet Rate changes are not counted.
Galtier N. and Jean-Marie A., Markov-modulated Markov chains and the covarion process of molecular evolution (2004). Journal of Computational Biology, 11:727-33.
Definition at line 46 of file MarkovModulatedSubstitutionModel.h.
|
inline |
Build a new MarkovModulatedSubstitutionModel object.
model | The substitution model to use. Can be of any alphabet type, and will be owned by this instance. |
nbRates | The number of rate classes |
normalizeRateChanges | Tells if the branch lengths must be computed in terms of rate and state NB: In most cases, this parameter should be set to false. |
prefix | The parameter namespace to be forwarded to the AbstractParametrizable constructor. changes instead of state change only. |
Definition at line 139 of file MarkovModulatedSubstitutionModel.h.
References bpp::AbstractParameterAliasable::addParameters_(), model_, and nestedPrefix_.
MarkovModulatedSubstitutionModel::MarkovModulatedSubstitutionModel | ( | const MarkovModulatedSubstitutionModel & | model | ) |
Definition at line 16 of file MarkovModulatedSubstitutionModel.cpp.
|
inlinevirtual |
Definition at line 172 of file MarkovModulatedSubstitutionModel.h.
|
pure virtualinherited |
Implemented in bpp::TS98, bpp::TransitionFromTransitionModel, bpp::RegisterRatesSubstitutionModel, bpp::MultinomialFromTransitionModel, bpp::InMixedSubstitutionModel, bpp::G2001, bpp::FromMixtureSubstitutionModel, bpp::AbstractTransitionModel, bpp::AbstractFromSubstitutionModelTransitionModel, and bpp::AbstractBiblioTransitionModel.
Referenced by bpp::AbstractBiblioTransitionModel::addRateParameter(), bpp::AbstractFromSubstitutionModelTransitionModel::addRateParameter(), bpp::FromMixtureSubstitutionModel::addRateParameter(), bpp::InMixedSubstitutionModel::addRateParameter(), bpp::MultinomialFromTransitionModel::addRateParameter(), bpp::TransitionFromTransitionModel::addRateParameter(), and bpp::BppOSubstitutionModelFormat::updateParameters_().
|
inlineoverridevirtual |
Implements bpp::BranchModelInterface.
Definition at line 177 of file MarkovModulatedSubstitutionModel.h.
References model_.
|
overridepure virtual |
Implements bpp::AbstractParameterAliasable.
Implemented in bpp::TS98, and bpp::G2001.
|
inlineoverridevirtual |
Implements bpp::TransitionModelInterface.
Definition at line 299 of file MarkovModulatedSubstitutionModel.h.
References compFreq_.
|
inlineoverridevirtual |
Implements bpp::TransitionModelInterface.
Definition at line 301 of file MarkovModulatedSubstitutionModel.h.
References compFreq_.
|
pure virtualinherited |
Implements bpp::BranchModelInterface.
Implemented in bpp::AbstractLkTransitionModel.
|
inlineoverridevirtualinherited |
Implements bpp::TransitionModelInterface.
Definition at line 64 of file AbstractSubstitutionModel.h.
References bpp::TransitionModelInterface::getd2Pij_dt2(), and bpp::AbstractLkTransitionModel::lik_.
|
inlineoverridevirtual |
Implements bpp::TransitionModelInterface.
Definition at line 222 of file MarkovModulatedSubstitutionModel.h.
References getd2Pij_dt2().
|
pure virtualinherited |
Implements bpp::BranchModelInterface.
Implemented in bpp::AbstractLkTransitionModel.
|
inlineoverridevirtualinherited |
Implements bpp::TransitionModelInterface.
Definition at line 47 of file AbstractSubstitutionModel.h.
References bpp::TransitionModelInterface::getdPij_dt(), and bpp::AbstractLkTransitionModel::lik_.
|
inlineoverridevirtual |
Implements bpp::TransitionModelInterface.
Definition at line 221 of file MarkovModulatedSubstitutionModel.h.
References getdPij_dt().
|
inlineoverridevirtual |
Tell if eigenValues and Vectors must be computed.
Implements bpp::SubstitutionModelInterface.
Definition at line 297 of file MarkovModulatedSubstitutionModel.h.
References eigenDecompose_.
|
inlineoverridevirtual |
Set if eigenValues and Vectors must be computed.
Implements bpp::SubstitutionModelInterface.
Definition at line 295 of file MarkovModulatedSubstitutionModel.h.
References eigenDecompose_.
|
inlineoverridevirtual |
Implements bpp::SubstitutionModelInterface.
Definition at line 199 of file MarkovModulatedSubstitutionModel.h.
References exchangeability_.
|
inlineoverridevirtual |
Tells the model that a parameter value has changed.
This updates the matrices consequently.
Reimplemented from bpp::AbstractParameterAliasable.
Reimplemented in bpp::G2001.
Definition at line 309 of file MarkovModulatedSubstitutionModel.h.
References model_, updateMatrices_(), and updateRatesModel_().
Referenced by bpp::G2001::fireParameterChanged().
|
inlineoverridevirtual |
Implements bpp::TransitionModelInterface.
Definition at line 216 of file MarkovModulatedSubstitutionModel.h.
References freq_.
|
inlineoverridevirtual |
Exception | if no FrequenceSet is associated to this model. |
Implements bpp::BranchModelInterface.
Definition at line 238 of file MarkovModulatedSubstitutionModel.h.
|
inlineoverridevirtual |
See Kosiol and Goldman (2005), Molecular Biology And Evolution 22(2) 193-9.
Implements bpp::SubstitutionModelInterface.
Definition at line 201 of file MarkovModulatedSubstitutionModel.h.
References generator_.
|
inlineoverridevirtual |
Implements bpp::BranchModelInterface.
Definition at line 179 of file MarkovModulatedSubstitutionModel.h.
References model_.
Referenced by getInitValue().
|
inlineoverridevirtual |
index | The model state. |
Implements bpp::BranchModelInterface.
Definition at line 189 of file MarkovModulatedSubstitutionModel.h.
References stateMap_.
|
inlineoverridevirtual |
index | The model state. |
Implements bpp::BranchModelInterface.
Definition at line 191 of file MarkovModulatedSubstitutionModel.h.
References stateMap_.
Referenced by getInitValue().
|
inlineoverridevirtual |
Implements bpp::BranchModelInterface.
Definition at line 187 of file MarkovModulatedSubstitutionModel.h.
References stateMap_.
|
inlineoverridevirtual |
Implements bpp::SubstitutionModelInterface.
Definition at line 214 of file MarkovModulatedSubstitutionModel.h.
References rightEigenVectors_.
|
overridevirtual |
Implements bpp::TransitionModelInterface.
Definition at line 177 of file MarkovModulatedSubstitutionModel.cpp.
References d2pijt_, eigenValues_, bpp::VectorTools::exp(), leftEigenVectors_, bpp::MatrixTools::mult(), rightEigenVectors_, and bpp::VectorTools::sqr().
Referenced by d2Pij_dt2().
|
overridevirtual |
Implements bpp::TransitionModelInterface.
Definition at line 171 of file MarkovModulatedSubstitutionModel.cpp.
References dpijt_, eigenValues_, bpp::VectorTools::exp(), leftEigenVectors_, bpp::MatrixTools::mult(), and rightEigenVectors_.
Referenced by dPij_dt().
|
inlineoverridevirtual |
Implements bpp::SubstitutionModelInterface.
Definition at line 207 of file MarkovModulatedSubstitutionModel.h.
References eigenValues_.
|
inlineoverridevirtual |
Implements bpp::TransitionModelInterface.
Definition at line 197 of file MarkovModulatedSubstitutionModel.h.
References freq_.
|
inlineoverrideprotectedvirtual |
Implements bpp::TransitionModelInterface.
Definition at line 329 of file MarkovModulatedSubstitutionModel.h.
References freq_.
|
inlineoverridevirtual |
Implements bpp::SubstitutionModelInterface.
Definition at line 208 of file MarkovModulatedSubstitutionModel.h.
References iEigenValues_.
|
overridevirtual |
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.
i | the index of the state in the model. |
state | An observed state in the sequence/site. |
IndexOutOfBoundsException | if array position is out of range. |
BadIntException | if states are not allowed in the associated alphabet. |
Implements bpp::BranchModelInterface.
Definition at line 185 of file MarkovModulatedSubstitutionModel.cpp.
References getAlphabet(), getAlphabetStateAsInt(), model_, nbRates_, and nbStates_.
|
inlineoverridevirtual |
Get the state in the model corresponding to a particular state in the alphabet.
code | The alphabet state to check. |
Implements bpp::BranchModelInterface.
Definition at line 195 of file MarkovModulatedSubstitutionModel.h.
References stateMap_.
|
inlineoverridevirtual |
Get the state in the model corresponding to a particular state in the alphabet.
code | The alphabet state to check. |
Implements bpp::BranchModelInterface.
Definition at line 193 of file MarkovModulatedSubstitutionModel.h.
References stateMap_.
|
pure virtualinherited |
Get the name of the model.
Implemented in bpp::WordSubstitutionModel, bpp::TS98, bpp::TransitionFromTransitionModel, bpp::RegisterRatesSubstitutionModel, bpp::RE08, bpp::WAG01, bpp::UserProteinSubstitutionModel, bpp::LLG08_UL3, bpp::LLG08_UL3::EmbeddedModel, bpp::LLG08_UL2, bpp::LLG08_UL2::EmbeddedModel, bpp::LLG08_EX3, bpp::LLG08_EX3::EmbeddedModel, bpp::LLG08_EX2, bpp::LLG08_EX2::EmbeddedModel, bpp::LLG08_EHO, bpp::LLG08_EHO::EmbeddedModel, bpp::LGL08_CAT, bpp::LGL08_CAT::EmbeddedModel, bpp::LG08, bpp::JTT92, bpp::JCprot, bpp::DSO78, bpp::Coala, bpp::POMO, bpp::OneChangeTransitionModel, bpp::OneChangeRegisterTransitionModel, bpp::YpR_Gen, bpp::YpR_Sym, bpp::TN93, bpp::T92, bpp::SSR, bpp::RN95s, bpp::RN95, bpp::L95, bpp::K80, bpp::JCnuc, bpp::GTR, bpp::gBGC, bpp::F81, bpp::MultinomialFromTransitionModel, bpp::MixtureOfTransitionModels, bpp::MixtureOfATransitionModel, bpp::KroneckerWordSubstitutionModel, bpp::InMixedSubstitutionModel, bpp::G2001, bpp::FromMixtureSubstitutionModel, bpp::EquiprobableSubstitutionModel, bpp::D1WalkSubstitutionModel, bpp::YNGP_M9, bpp::YNGP_M8, bpp::YNGP_M7, bpp::YNGP_M3, bpp::YNGP_M2, bpp::YNGP_M10, bpp::YNGP_M1, bpp::YN98, bpp::TripletSubstitutionModel, bpp::SENCA, bpp::RELAX, bpp::MG94, bpp::KroneckerCodonDistanceSubstitutionModel, bpp::KroneckerCodonDistanceFrequenciesSubstitutionModel, bpp::KCM, bpp::GY94, bpp::DFPDistanceFrequenciesSubstitutionModel, bpp::DFP07, bpp::CodonSameAARateSubstitutionModel, bpp::CodonDistancePhaseFrequenciesSubstitutionModel, bpp::CodonDistanceFrequenciesSubstitutionModel, bpp::CodonAdHocSubstitutionModel, bpp::AnonymousSubstitutionModel, bpp::AbstractWrappedModel, bpp::TwoParameterBinarySubstitutionModel, bpp::LG10_EX_EHO, bpp::LG10_EX_EHO::EmbeddedModel, bpp::HKY85, bpp::F84, bpp::CodonDistanceSubstitutionModel, and bpp::BinarySubstitutionModel.
Referenced by bpp::YpR::checkModel(), bpp::TransitionMatrixFromModel::compute(), bpp::AbstractWrappedModel::getName(), bpp::FromMixtureSubstitutionModel::getName(), bpp::AbstractSubstitutionModel::getPij_t(), bpp::BppOSubstitutionModelFormat::initialize_(), bpp::MixtureOfSubstitutionModels::MixtureOfSubstitutionModels(), bpp::MixtureOfATransitionModel::model(), bpp::LegacyPhylogeneticsApplicationTools::setSubstitutionModelParametersInitialValuesWithAliases(), bpp::ModelPath::toString(), and bpp::BppOSubstitutionModelFormat::write().
|
inlineoverridevirtual |
Get the number of states.
For most models, this equals the size of the alphabet.
Implements bpp::BranchModelInterface.
Definition at line 181 of file MarkovModulatedSubstitutionModel.h.
References stateMap_.
Referenced by setDiagonal().
|
overridevirtual |
Implements bpp::TransitionModelInterface.
Definition at line 162 of file MarkovModulatedSubstitutionModel.cpp.
References eigenValues_, bpp::VectorTools::exp(), leftEigenVectors_, bpp::MatrixTools::mult(), nbRates_, nbStates_, pijt_, and rightEigenVectors_.
Referenced by Pij_t().
|
inlineoverridevirtual |
Get the rate.
Implements bpp::BranchModelInterface.
Reimplemented in bpp::TS98.
Definition at line 260 of file MarkovModulatedSubstitutionModel.h.
References model_.
|
inline |
Get the rate category corresponding to a particular state in the compound model.
i | The state. |
Definition at line 255 of file MarkovModulatedSubstitutionModel.h.
References nbStates_.
|
inlineoverridevirtual |
Implements bpp::SubstitutionModelInterface.
Definition at line 213 of file MarkovModulatedSubstitutionModel.h.
References leftEigenVectors_.
|
inlineoverridevirtual |
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.
Implements bpp::SubstitutionModelInterface.
Definition at line 282 of file MarkovModulatedSubstitutionModel.h.
References bpp::MatrixTools::diag(), freq_, and generator_.
Referenced by updateMatrices_().
|
inlineoverridevirtual |
Implements bpp::BranchModelInterface.
Definition at line 185 of file MarkovModulatedSubstitutionModel.h.
References stateMap_.
|
inlineoverridevirtual |
Implements bpp::SubstitutionModelInterface.
Definition at line 210 of file MarkovModulatedSubstitutionModel.h.
|
inlineoverridevirtual |
Implements bpp::SubstitutionModelInterface.
Definition at line 211 of file MarkovModulatedSubstitutionModel.h.
|
inlineoverridevirtual |
returns if model is scalable
Implements bpp::SubstitutionModelInterface.
Definition at line 264 of file MarkovModulatedSubstitutionModel.h.
References model_.
Referenced by updateMatrices_().
|
pure virtualinherited |
This method is used to compute likelihoods in recursions. It computes the probability of a vector given a start state.
values | An vector of states on the site. |
t | the branch length |
Implements bpp::BranchModelInterface.
Implemented in bpp::AbstractLkTransitionModel.
|
inlineoverridevirtualinherited |
This method is used to compute likelihoods in recursions. It computes the probability of a vector given a start state.
values | An vector of states on the site. |
t | the branch length |
Implements bpp::TransitionModelInterface.
Definition at line 30 of file AbstractSubstitutionModel.h.
References bpp::TransitionModelInterface::getPij_t(), and bpp::AbstractLkTransitionModel::lik_.
|
inline |
Definition at line 243 of file MarkovModulatedSubstitutionModel.h.
References model_.
Referenced by bpp::BppOSubstitutionModelFormat::write().
|
inlineoverridevirtual |
Normalize the generator.
Implements bpp::SubstitutionModelInterface.
Definition at line 274 of file MarkovModulatedSubstitutionModel.h.
References model_, and updateMatrices_().
MarkovModulatedSubstitutionModel & MarkovModulatedSubstitutionModel::operator= | ( | const MarkovModulatedSubstitutionModel & | model | ) |
Definition at line 43 of file MarkovModulatedSubstitutionModel.cpp.
References compFreq_, d2pijt_, dpijt_, eigenDecompose_, eigenValues_, exchangeability_, freq_, generator_, iEigenValues_, leftEigenVectors_, model_, nbRates_, nbStates_, nestedPrefix_, normalizeRateChanges_, pijt_, rates_, ratesExchangeability_, ratesFreq_, ratesGenerator_, rightEigenVectors_, and stateMap_.
Referenced by bpp::G2001::operator=().
|
inlineoverridevirtual |
Implements bpp::TransitionModelInterface.
Definition at line 220 of file MarkovModulatedSubstitutionModel.h.
References getPij_t().
|
inlineoverridevirtual |
A method for computing all necessary matrices.
Implements bpp::SubstitutionModelInterface.
Definition at line 218 of file MarkovModulatedSubstitutionModel.h.
References generator_.
|
overridevirtual |
set the diagonal of the generator such that sum on each line equals 0.
Implements bpp::SubstitutionModelInterface.
Definition at line 144 of file MarkovModulatedSubstitutionModel.cpp.
References generator_, getNumberOfStates(), and RowMatrix< double >::getRow().
|
inlineoverridevirtual |
Set equilibrium frequencies.
frequencies | The map of the frequencies to use. |
Implements bpp::TransitionModelInterface.
Definition at line 232 of file MarkovModulatedSubstitutionModel.h.
References model_, and updateMatrices_().
|
inlineoverridevirtual |
Set equilibrium frequencies equal to the frequencies estimated from the data.
data | The sequences to use. |
pseudoCount | A quantity 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
|
Implements bpp::TransitionModelInterface.
Definition at line 226 of file MarkovModulatedSubstitutionModel.h.
References model_, and updateMatrices_().
|
overridevirtual |
Reimplemented from bpp::AbstractParameterAliasable.
Definition at line 202 of file MarkovModulatedSubstitutionModel.cpp.
References model_, nestedPrefix_, and bpp::AbstractParameterAliasable::setNamespace().
Referenced by bpp::G2001::setNamespace().
|
inlineoverridevirtual |
Set the rate of the model (must be positive).
rate | must be positive. |
Implements bpp::BranchModelInterface.
Reimplemented in bpp::TS98.
Definition at line 262 of file MarkovModulatedSubstitutionModel.h.
References model_.
|
inlineoverridevirtual |
sets if model is scalable, ie scale can be changed. Default : true, set to false to avoid normalization for example.
Implements bpp::SubstitutionModelInterface.
Definition at line 269 of file MarkovModulatedSubstitutionModel.h.
References model_.
|
inlineoverridevirtual |
Multiplies the current generator by the given scale.
scale | the scale by which the generator is multiplied. |
Implements bpp::SubstitutionModelInterface.
Definition at line 289 of file MarkovModulatedSubstitutionModel.h.
References model_, and updateMatrices_().
Referenced by updateMatrices_().
|
inlineoverridevirtual |
By definition Sij(i,j) = Sij(j,i).
Implements bpp::SubstitutionModelInterface.
Definition at line 217 of file MarkovModulatedSubstitutionModel.h.
References exchangeability_.
|
inlineoverridevirtual |
Implements bpp::BranchModelInterface.
Definition at line 183 of file MarkovModulatedSubstitutionModel.h.
References stateMap_.
|
protectedvirtual |
Definition at line 74 of file MarkovModulatedSubstitutionModel.cpp.
References bpp::MatrixTools::add(), d2pijt_, bpp::MatrixTools::diag(), dpijt_, eigenValues_, exchangeability_, freq_, generator_, RowMatrix< double >::getNumberOfColumns(), bpp::EigenValue< class >::getRealEigenValues(), getScale(), bpp::EigenValue< class >::getV(), iEigenValues_, bpp::MatrixTools::inv(), isScalable(), bpp::MatrixTools::kroneckerMult(), bpp::VectorTools::kroneckerMult(), leftEigenVectors_, model_, bpp::MatrixTools::mult(), nbRates_, nbStates_, normalizeRateChanges_, pijt_, rates_, ratesExchangeability_, ratesFreq_, ratesGenerator_, RowMatrix< double >::resize(), rightEigenVectors_, bpp::MatrixTools::scale(), and setScale().
Referenced by fireParameterChanged(), bpp::G2001::G2001(), normalize(), setFreq(), setFreqFromData(), setScale(), and bpp::TS98::TS98().
|
protectedpure virtual |
Update the rates vector, generator and equilibrium frequencies.
This method must be implemented by the derived class. It is called by the fireParameterChanged() method.
Implemented in bpp::TS98, and bpp::G2001.
Referenced by fireParameterChanged().
|
protected |
Tell if the equilibrium frequencies should be computed from the generator.
Definition at line 110 of file MarkovModulatedSubstitutionModel.h.
Referenced by computeFrequencies(), and operator=().
|
mutableprotected |
Definition at line 117 of file MarkovModulatedSubstitutionModel.h.
Referenced by getd2Pij_dt2(), operator=(), and updateMatrices_().
|
mutableprotected |
Definition at line 116 of file MarkovModulatedSubstitutionModel.h.
Referenced by getdPij_dt(), operator=(), and updateMatrices_().
|
protected |
Tell if the eigen decomposition should be performed.
Definition at line 103 of file MarkovModulatedSubstitutionModel.h.
Referenced by enableEigenDecomposition(), and operator=().
|
protected |
The vector of real parts of eigen values.
Definition at line 92 of file MarkovModulatedSubstitutionModel.h.
Referenced by getd2Pij_dt2(), getdPij_dt(), getEigenValues(), getPij_t(), operator=(), and updateMatrices_().
|
protected |
The exchangeability matrix of the model.
Definition at line 77 of file MarkovModulatedSubstitutionModel.h.
Referenced by exchangeabilityMatrix(), operator=(), Sij(), and updateMatrices_().
|
protected |
The vector of equilibrium frequencies.
Definition at line 122 of file MarkovModulatedSubstitutionModel.h.
Referenced by freq(), getFrequencies(), getFrequencies_(), getScale(), operator=(), and updateMatrices_().
|
protected |
The generator matrix of the model.
Definition at line 72 of file MarkovModulatedSubstitutionModel.h.
Referenced by generator(), getScale(), operator=(), Qij(), setDiagonal(), and updateMatrices_().
|
protected |
The vector of imaginary parts of the eigen values (zero in case of reversible pmodel).
Definition at line 98 of file MarkovModulatedSubstitutionModel.h.
Referenced by getIEigenValues(), operator=(), and updateMatrices_().
|
protected |
The matrix made of left eigen vectors (by row).
Definition at line 82 of file MarkovModulatedSubstitutionModel.h.
Referenced by getd2Pij_dt2(), getdPij_dt(), getPij_t(), getRowLeftEigenVectors(), operator=(), and updateMatrices_().
|
mutableprivateinherited |
Definition at line 23 of file AbstractSubstitutionModel.h.
Referenced by bpp::AbstractLkTransitionModel::d2Lik_dt2(), bpp::AbstractLkTransitionModel::dLik_dt(), and bpp::AbstractLkTransitionModel::Lik_t().
|
protected |
Definition at line 52 of file MarkovModulatedSubstitutionModel.h.
Referenced by alphabet(), fireParameterChanged(), getAlphabet(), getInitValue(), getRate(), isScalable(), MarkovModulatedSubstitutionModel(), nestedModel(), normalize(), operator=(), setFreq(), setFreqFromData(), setNamespace(), setRate(), setScalable(), setScale(), and updateMatrices_().
|
protected |
Definition at line 55 of file MarkovModulatedSubstitutionModel.h.
Referenced by bpp::G2001::G2001(), getInitValue(), getPij_t(), operator=(), updateMatrices_(), and bpp::G2001::updateRatesModel_().
|
protected |
Definition at line 54 of file MarkovModulatedSubstitutionModel.h.
Referenced by getInitValue(), getPij_t(), getRate(), operator=(), and updateMatrices_().
|
protected |
Definition at line 126 of file MarkovModulatedSubstitutionModel.h.
Referenced by MarkovModulatedSubstitutionModel(), operator=(), and setNamespace().
|
protected |
Definition at line 124 of file MarkovModulatedSubstitutionModel.h.
Referenced by operator=(), and updateMatrices_().
|
mutableprotected |
These ones are for bookkeeping:
Definition at line 115 of file MarkovModulatedSubstitutionModel.h.
Referenced by getPij_t(), operator=(), and updateMatrices_().
|
protected |
Definition at line 63 of file MarkovModulatedSubstitutionModel.h.
Referenced by operator=(), updateMatrices_(), bpp::G2001::updateRatesModel_(), and bpp::TS98::updateRatesModel_().
|
protected |
Definition at line 64 of file MarkovModulatedSubstitutionModel.h.
Referenced by operator=(), updateMatrices_(), bpp::G2001::updateRatesModel_(), and bpp::TS98::updateRatesModel_().
|
protected |
Definition at line 65 of file MarkovModulatedSubstitutionModel.h.
Referenced by bpp::G2001::G2001(), operator=(), updateMatrices_(), and bpp::TS98::updateRatesModel_().
|
protected |
Definition at line 67 of file MarkovModulatedSubstitutionModel.h.
Referenced by operator=(), and updateMatrices_().
|
protected |
The matrix made of right eigen vectors (by column).
Definition at line 87 of file MarkovModulatedSubstitutionModel.h.
Referenced by getColumnRightEigenVectors(), getd2Pij_dt2(), getdPij_dt(), getPij_t(), operator=(), and updateMatrices_().
|
protected |
Definition at line 53 of file MarkovModulatedSubstitutionModel.h.
Referenced by getAlphabetStateAsChar(), getAlphabetStateAsInt(), getAlphabetStates(), getModelStates(), getNumberOfStates(), getStateMap(), operator=(), and stateMap().