5#ifndef BPP_PHYL_MODEL_ABSTRACTSUBSTITUTIONMODEL_H
6#define BPP_PHYL_MODEL_ABSTRACTSUBSTITUTIONMODEL_H
23 mutable Eigen::VectorXd
lik_;
30 const Eigen::VectorXd&
Lik_t(
const Eigen::VectorXd& values,
double t)
const override
32 lik_ = Eigen::VectorXd::Zero(values.size());
36 for (
auto i = 0; i < values.size(); ++i)
38 for (
auto j = 0; j < values.size(); ++j)
40 lik_(i) += pij(
size_t(i),
size_t(j)) * values(j);
47 const Eigen::VectorXd&
dLik_dt(
const Eigen::VectorXd& values,
double t)
const override
49 lik_ = Eigen::VectorXd::Zero(values.size());
53 for (
auto i = 0; i < values.size(); ++i)
55 for (
auto j = 0; j < values.size(); ++j)
57 lik_(i) += pij(
size_t(i),
size_t(j)) * values(j);
64 const Eigen::VectorXd&
d2Lik_dt2(
const Eigen::VectorXd& values,
double t)
const override
66 lik_ = Eigen::VectorXd::Zero(values.size());
70 for (
auto i = 0; i < values.size(); ++i)
72 for (
auto j = 0; j < values.size(); ++j)
74 lik_(i) += pij(
size_t(i),
size_t(j)) * values(j);
156 std::shared_ptr<const Alphabet> alpha,
157 std::shared_ptr<const StateMapInterface>
stateMap,
158 const std::string& prefix);
195 virtual double freq(
size_t i)
const override {
return freq_[i]; }
197 virtual double Pij_t (
size_t i,
size_t j,
double t)
const override {
return getPij_t(t) (i, j); }
198 virtual double dPij_dt (
size_t i,
size_t j,
double t)
const override {
return getdPij_dt(t) (i, j); }
201 double getInitValue(
size_t i,
int state)
const override;
205 virtual void setFreq(std::map<int, double>& freqs)
override;
209 throw Exception(
"TransitionModel::frequencySet(). No associated FrequencySet object.");
225 if (parameters.
size() != 1)
278 virtual double getRate()
const override;
280 virtual void setRate(
double rate)
override;
361 std::shared_ptr<const Alphabet> alpha,
362 std::shared_ptr<const StateMapInterface>
stateMap,
363 const std::string& prefix);
530 std::shared_ptr<const Alphabet> alpha,
531 std::shared_ptr<const StateMapInterface>
stateMap,
532 const std::string& prefix) :
540 for (
auto& fr :
freq_)
542 fr = 1.0 /
static_cast<double>(
size_);
Virtual class of a Transition Model related to a given SubstitutionModel.
Partial implementation of the TransitionModel interface, with function for likelihood computations.
const Eigen::VectorXd & dLik_dt(const Eigen::VectorXd &values, double t) const override
virtual ~AbstractLkTransitionModel()
AbstractLkTransitionModel()
const Eigen::VectorXd & Lik_t(const Eigen::VectorXd &values, double t) const override
const Eigen::VectorXd & d2Lik_dt2(const Eigen::VectorXd &values, double t) const override
virtual void fireParameterChanged(const ParameterList ¶meters)
std::string getNamespace() const override
Partial implementation of the ReversibleSubstitutionModel interface.
AbstractReversibleSubstitutionModel(std::shared_ptr< const Alphabet > alpha, std::shared_ptr< const StateMapInterface > stateMap, const std::string &prefix)
virtual ~AbstractReversibleSubstitutionModel()
virtual AbstractReversibleSubstitutionModel * clone() const override=0
virtual void updateMatrices_() override
Compute and diagonalize the matrix, and fill the eigenValues_, leftEigenVectors_ and rightEigenVecto...
bool computeFreq_
if the Frequencies must be computed from the generator
RowMatrix< double > generator_
The generator matrix of the model.
const Matrix< double > & generator() const
bool isDiagonalizable_
boolean value for diagonalizability in R of the generator_
bool isDiagonalizable() const
bool isNonSingular_
boolean value for non-singularity of rightEigenVectors_
std::vector< RowMatrix< double > > vPowGen_
vector of the powers of generator_ for Taylor development (if rightEigenVectors_ is singular).
Vdouble eigenValues_
The vector of eigen values.
AbstractSubstitutionModel(std::shared_ptr< const Alphabet > alpha, std::shared_ptr< const StateMapInterface > stateMap, const std::string &prefix)
AbstractSubstitutionModel & operator=(const AbstractSubstitutionModel &model)
void normalize()
normalize the generator
virtual ~AbstractSubstitutionModel()
void setDiagonal()
set the diagonal of the generator such that sum on each line equals 0.
bool eigenDecompose_
Tell if the eigen decomposition should be performed.
const Matrix< double > & exchangeabilityMatrix() const
const Vdouble & getEigenValues() const
virtual double Qij(size_t i, size_t j) const
A method for computing all necessary matrices.
void setScalable(bool scalable)
sets if model is scalable, ie scale can be changed. Default : true, set to false to avoid normalizati...
bool computeFrequencies() const
void enableEigenDecomposition(bool yn)
Set if eigenValues and Vectors must be computed.
void computeFrequencies(bool yn)
Set if equilibrium frequencies should be computed.
RowMatrix< double > exchangeability_
The exchangeability matrix of the model, defined as . When the model is reversible,...
const Matrix< double > & getdPij_dt(double t) const
Vdouble iEigenValues_
The vector of the imaginary part of the eigen values.
RowMatrix< double > leftEigenVectors_
The matrix made of left eigen vectors (by row) if rightEigenVectors_ is non-singular.
bool isScalable_
If the model is scalable (ie generator can be normalized automatically).
AbstractSubstitutionModel(const AbstractSubstitutionModel &model)
const Matrix< double > & getPij_t(double t) const
virtual void updateMatrices_()
Diagonalize the matrix, and fill the eigenValues_, iEigenValues_, leftEigenVectors_ and rightEigenVe...
bool isNonSingular() const
double Sij(size_t i, size_t j) const
virtual bool isScalable() const
returns if model is scalable
const Vdouble & getIEigenValues() const
RowMatrix< double > rightEigenVectors_
The matrix made of right eigen vectors (by column).
bool enableEigenDecomposition()
Tell if eigenValues and Vectors must be computed.
const Matrix< double > & getRowLeftEigenVectors() const
const Matrix< double > & getColumnRightEigenVectors() const
const Matrix< double > & getd2Pij_dt2(double t) const
void setScale(double scale)
Multiplies the current generator by the given scale.
double getScale() const
return scale
RowMatrix< double > tmpMat_
For computational issues.
Partial implementation of the TransitionModel interface.
void setFreqFromData(const SequenceDataInterface &data, double pseudoCount=0) override
Set equilibrium frequencies equal to the frequencies estimated from the data.
void addRateParameter() override
add a "rate" parameter to the model, that handles the overall rate of the process.
std::string getAlphabetStateAsChar(size_t index) const override
RowMatrix< double > d2pijt_
std::shared_ptr< const StateMapInterface > stateMap_
The map of model states with alphabet states.
const Alphabet & alphabet() const override
virtual ~AbstractTransitionModel()
virtual double d2Pij_dt2(size_t i, size_t j, double t) const override
AbstractTransitionModel(const AbstractTransitionModel &model)=default
virtual void fireParameterChanged(const ParameterList ¶meters) override
Tells the model that a parameter value has changed.
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.
const Vdouble & getFrequencies() const override
RowMatrix< double > pijt_
These ones are for bookkeeping:
size_t size_
The number of states.
Vdouble freq_
The vector of equilibrium frequencies.
virtual void setRate(double rate) override
Set the rate of the model (must be positive).
virtual const Matrix< double > & getdPij_dt(double t) const override=0
std::shared_ptr< const Alphabet > getAlphabet() const override
AbstractTransitionModel & operator=(const AbstractTransitionModel &model)=default
int getAlphabetStateAsInt(size_t index) const override
virtual const Matrix< double > & getPij_t(double t) const override=0
virtual void updateMatrices_()=0
Diagonalize the matrix, and fill the eigenValues_, iEigenValues_, leftEigenVectors_ and rightEigenVe...
RowMatrix< double > dpijt_
double rate_
The rate of the model (default: 1). The generator (and all its vectorial components) is independent o...
std::vector< size_t > getModelStates(int code) const override
Get the state in the model corresponding to a particular state in the alphabet.
virtual double freq(size_t i) const override
virtual const Matrix< double > & getd2Pij_dt2(double t) const override=0
std::shared_ptr< const StateMapInterface > getStateMap() const override
std::shared_ptr< const Alphabet > alphabet_
The alphabet relevant to this model.
void setVerboseLevel(short level)
const StateMapInterface & stateMap() const override
const std::vector< int > & getAlphabetStates() const override
short verboseLevel() const
size_t getNumberOfStates() const override
Get the number of states.
Vdouble & getFrequencies_() override
virtual double getRate() const override
The rate of the substitution process.
double getInitValue(size_t i, int state) const override
bool computeFrequencies() const override
AbstractTransitionModel(std::shared_ptr< const Alphabet > alpha, std::shared_ptr< const StateMapInterface > stateMap, const std::string &prefix)
virtual double Pij_t(size_t i, size_t j, double t) const override
virtual void setFreq(std::map< int, double > &freqs) override
Set equilibrium frequencies.
const FrequencySetInterface & frequencySet() const override
virtual double dPij_dt(size_t i, size_t j, double t) const override
Parametrize a set of state frequencies.
SubModel taken from a MixedTransitionModel, kept in the context of the MixedTransitionModel (see From...
From a model, compute transition probabilities given there is at least a change of a category (ie a n...
virtual bool hasParameter(const std::string &name) const
virtual double getParameterValue(const std::string &name) const
Interface for reversible substitution models.
Map the states of a given alphabet which have a model state.
Interface for all substitution models.
Interface for all transition models.
virtual const Matrix< double > & getPij_t(double t) const =0
virtual const Matrix< double > & getdPij_dt(double t) const =0
virtual const Matrix< double > & getd2Pij_dt2(double t) const =0
Defines the basic types of data flow nodes.
std::vector< double > Vdouble