bpp-phyl3  3.0.0
HmmLikelihood_DF.h
Go to the documentation of this file.
1 // SPDX-FileCopyrightText: The Bio++ Development Group
2 //
3 // SPDX-License-Identifier: CECILL-2.1
4 
5 #ifndef BPP_PHYL_LIKELIHOOD_PHYLOLIKELIHOODS_HMMLIKELIHOOD_DF_H
6 #define BPP_PHYL_LIKELIHOOD_PHYLOLIKELIHOODS_HMMLIKELIHOOD_DF_H
7 
10 #include <Bpp/Numeric/NumTools.h>
11 
12 #include "../DataFlow/LikelihoodCalculation.h"
13 #include "../DataFlow/TransitionMatrix.h"
16 
17 // From the STL:
18 #include <vector>
19 #include <memory>
20 
21 namespace bpp
22 {
34 {
35 private:
37 
38 protected:
43  std::shared_ptr<HmmStateAlphabet> hiddenAlphabet_;
44  std::shared_ptr<HmmPhyloEmissionProbabilities> emissionProbabilities_;
45 
46  /************************
47  * DF objects & their targets
48  *
49  */
50 
56  std::shared_ptr<ConfiguredTransitionMatrix> matrix_;
57 
64 
71 
78 
85 
86 
96 
106 
107  Eigen::Index nbStates_, nbSites_;
108 
109 public:
120  Context& context,
121  std::shared_ptr<HmmStateAlphabet> hiddenAlphabet,
122  std::shared_ptr<HmmTransitionMatrix> transitionMatrix,
123  std::shared_ptr<HmmPhyloEmissionProbabilities> emissionProbabilities,
124  const std::string& prefix = "");
125 
128  context_(lik.context_),
130  matrix_(lik.matrix_),
131  hmmEq_(lik.hmmEq_),
132  hmmTrans_(lik.hmmTrans_),
133  hmmEmis_(lik.hmmEmis_),
137  nbStates_(lik.nbStates_),
138  nbSites_(lik.nbSites_)
139  {}
140 
141  virtual ~HmmLikelihood_DF() {}
142 
143  HmmLikelihood_DF* clone() const { return new HmmLikelihood_DF(*this); }
144 
145  void makeLikelihoods() {}
146 
147 public:
149 
151 
152  /*
153  *@ brief Access to the Transition Matrix
154  *
155  * !! No check if DF up to date
156  *
157  */
158  const Eigen::MatrixXd& getHmmTransitionMatrix() const { return hmmTrans_->accessValueConst(); }
159 
161 
163 
164  void setParameters(const ParameterList& pl)
165  {
167  }
168 
169  void setNamespace(const std::string& nameSpace);
170 
171  const Eigen::MatrixXd& getHiddenStatesPosteriorProbabilities() const
172  {
173  return hiddenPostProb_->targetValue();
174  }
175 
176  Eigen::VectorXd getHiddenStatesPosteriorProbabilitiesForASite(size_t site) const
177  {
178  auto& mat = hiddenPostProb_->targetValue();
179  return mat.col(Eigen::Index(site));
180  }
181 
182  DataLik getLikelihoodForASite(size_t site) const
183  {
185 
186  return hmmEmis_->accessValueConst().col(int(site)).dot(vec);
187  }
188 };
189 }
190 #endif // BPP_PHYL_LIKELIHOOD_PHYLOLIKELIHOODS_HMMLIKELIHOOD_DF_H
void setParametersValues(const ParameterList &parameters) override
Context for dataflow node construction.
Definition: DataFlow.h:527
A simple implementation of hidden Markov models recursion, in DataFlow implementation.
Eigen::VectorXd getHiddenStatesPosteriorProbabilitiesForASite(size_t site) const
const Eigen::MatrixXd & getHiddenStatesPosteriorProbabilities() const
ValueRef< Eigen::MatrixXd > hmmTrans_
void setNamespace(const std::string &nameSpace)
HmmLikelihood_DF * clone() const
HmmLikelihood_DF(Context &context, std::shared_ptr< HmmStateAlphabet > hiddenAlphabet, std::shared_ptr< HmmTransitionMatrix > transitionMatrix, std::shared_ptr< HmmPhyloEmissionProbabilities > emissionProbabilities, const std::string &prefix="")
Build a new HmmLikelihood_DF object.
ValueRef< RowLik > forwardLik_
std::shared_ptr< HmmStateAlphabet > hiddenAlphabet_
The alphabet describing the hidden states.
DataLik getLikelihoodForASite(size_t site) const
const Eigen::MatrixXd & getHmmTransitionMatrix() const
std::shared_ptr< HmmPhyloEmissionProbabilities > emissionProbabilities_
HmmPhyloEmissionProbabilities & getHmmEmissionProbabilities()
ValueRef< MatrixLik > hmmEmis_
const HmmPhyloEmissionProbabilities & getHmmEmissionProbabilities() const
ValueRef< Eigen::MatrixXd > hiddenPostProb_
void setParameters(const ParameterList &pl)
const HmmStateAlphabet & hmmStateAlphabet() const
ValueRef< Eigen::MatrixXd > backwardLik_
std::shared_ptr< ConfiguredTransitionMatrix > matrix_
ValueRef< Eigen::VectorXd > hmmEq_
HmmLikelihood_DF(const HmmLikelihood_DF &lik)
HmmStateAlphabet & hmmStateAlphabet()
Emission probabilities in the context of DF phylolikeihoods.
Defines the basic types of data flow nodes.
std::shared_ptr< Value< T > > ValueRef
Shared pointer alias for Value<T>.
Definition: DataFlow.h:84