bpp-core3  3.0.0
RescaledHmmLikelihood.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_NUMERIC_HMM_RESCALEDHMMLIKELIHOOD_H
6 #define BPP_NUMERIC_HMM_RESCALEDHMMLIKELIHOOD_H
7 
8 
9 #include "../AbstractParametrizable.h"
10 #include "../Matrix/Matrix.h"
11 #include "HmmLikelihood.h"
12 
13 // From the STL:
14 #include <vector>
15 #include <memory>
16 
17 namespace bpp
18 {
26  public virtual AbstractHmmLikelihood,
28 {
29 private:
33  std::shared_ptr<HmmStateAlphabet> hiddenAlphabet_;
34  std::shared_ptr<HmmTransitionMatrix> transitionMatrix_;
35  std::shared_ptr<HmmEmissionProbabilities> emissionProbabilities_;
36 
49  std::vector<double> likelihood_;
50 
58  mutable std::vector<std::vector<double>> dLikelihood_;
59  mutable std::vector<std::vector<double>> d2Likelihood_;
60 
68  mutable std::vector<std::vector<double>> backLikelihood_;
70 
78  std::vector<double> scales_;
79 
80  mutable std::vector<double> dScales_;
81  mutable std::vector<double> d2Scales_;
82  double logLik_;
83 
84  std::vector<size_t> breakPoints_;
85 
87 
88 public:
98  std::shared_ptr<HmmStateAlphabet> hiddenAlphabet,
99  std::shared_ptr<HmmTransitionMatrix> transitionMatrix,
100  std::shared_ptr<HmmEmissionProbabilities> emissionProbabilities,
101  const std::string& prefix);
102 
106  hiddenAlphabet_(dynamic_cast<HmmStateAlphabet*>(lik.hiddenAlphabet_->clone())),
107  transitionMatrix_(dynamic_cast<HmmTransitionMatrix*>(lik.transitionMatrix_->clone())),
108  emissionProbabilities_(dynamic_cast<HmmEmissionProbabilities*>(lik.emissionProbabilities_->clone())),
109  likelihood_(lik.likelihood_),
110  dLikelihood_(lik.dLikelihood_),
111  d2Likelihood_(lik.d2Likelihood_),
112  backLikelihood_(lik.backLikelihood_),
113  backLikelihoodUpToDate_(lik.backLikelihoodUpToDate_),
114  scales_(lik.scales_),
115  dScales_(lik.dScales_),
116  d2Scales_(lik.d2Scales_),
117  logLik_(lik.logLik_),
118  breakPoints_(lik.breakPoints_),
119  nbStates_(lik.nbStates_),
120  nbSites_(lik.nbSites_)
121  {
122  // Now adjust pointers:
123  transitionMatrix_->setHmmStateAlphabet(hiddenAlphabet_);
124  emissionProbabilities_->setHmmStateAlphabet(hiddenAlphabet_);
125  }
126 
128  {
130  AbstractParametrizable::operator=(lik);
131  hiddenAlphabet_ = std::unique_ptr<HmmStateAlphabet>(dynamic_cast<HmmStateAlphabet*>(lik.hiddenAlphabet_->clone()));
132  transitionMatrix_ = std::unique_ptr<HmmTransitionMatrix>(dynamic_cast<HmmTransitionMatrix*>(lik.transitionMatrix_->clone()));
133  emissionProbabilities_ = std::unique_ptr<HmmEmissionProbabilities>(dynamic_cast<HmmEmissionProbabilities*>(lik.emissionProbabilities_->clone()));
134  likelihood_ = lik.likelihood_;
135  dLikelihood_ = lik.dLikelihood_;
136  d2Likelihood_ = lik.d2Likelihood_;
137  backLikelihood_ = lik.backLikelihood_;
138  backLikelihoodUpToDate_ = lik.backLikelihoodUpToDate_;
139  scales_ = lik.scales_;
140  dScales_ = lik.dScales_;
141  d2Scales_ = lik.d2Scales_;
142  logLik_ = lik.logLik_;
143  breakPoints_ = lik.breakPoints_;
144  nbStates_ = lik.nbStates_;
145  nbSites_ = lik.nbSites_;
146 
147  // Now adjust pointers:
148  transitionMatrix_->setHmmStateAlphabet(hiddenAlphabet_);
149  emissionProbabilities_->setHmmStateAlphabet(hiddenAlphabet_);
150  return *this;
151  }
152 
154 
155  RescaledHmmLikelihood* clone() const override { return new RescaledHmmLikelihood(*this); }
156 
157 public:
158  const HmmStateAlphabet& hmmStateAlphabet() const override { return *hiddenAlphabet_; }
159  std::shared_ptr<const HmmStateAlphabet> getHmmStateAlphabet() const override { return hiddenAlphabet_; }
160 
162  std::shared_ptr<HmmStateAlphabet> getHmmStateAlphabet() override { return hiddenAlphabet_; }
163 
164  const HmmTransitionMatrix& hmmTransitionMatrix() const override { return *transitionMatrix_; }
165  std::shared_ptr<const HmmTransitionMatrix> getHmmTransitionMatrix() const override { return transitionMatrix_; }
166 
168  std::shared_ptr<HmmTransitionMatrix> getHmmTransitionMatrix() override { return transitionMatrix_; }
169 
171  std::shared_ptr<const HmmEmissionProbabilities> getHmmEmissionProbabilities() const override { return emissionProbabilities_; }
172 
174  std::shared_ptr<HmmEmissionProbabilities> getHmmEmissionProbabilities() override { return emissionProbabilities_; }
175 
176  void setBreakPoints(const std::vector<size_t>& breakPoints) override
177  {
178  breakPoints_ = breakPoints;
179  computeForward_();
180  backLikelihoodUpToDate_ = false;
181  }
182 
183  const std::vector<size_t>& getBreakPoints() const override { return breakPoints_; }
184 
185  void setParameters(const ParameterList& pl) override
186  {
188  }
189 
190  double getValue() const override { return -logLik_; }
191 
192  double getLogLikelihood() const override { return logLik_; }
193 
194  double getLikelihoodForASite(size_t site) const override;
195 
196  double getDLogLikelihoodForASite(size_t site) const override;
197 
198  double getD2LogLikelihoodForASite(size_t site) const override;
199 
200  Vdouble getLikelihoodForEachSite() const override;
201 
202  void setNamespace(const std::string& nameSpace) override;
203 
204  void fireParameterChanged(const ParameterList& pl) override;
205 
206  Vdouble getHiddenStatesPosteriorProbabilitiesForASite(size_t site) const override;
207 
208  void getHiddenStatesPosteriorProbabilities(std::vector< std::vector<double>>& probs, bool append = false) const override;
209 
210 protected:
211  void computeForward_();
212  void computeBackward_() const;
213 
214 
215  void computeDLikelihood_() const override
216  {
218  }
219 
220  void computeD2Likelihood_() const override
221  {
223  }
224 
225  void computeDForward_() const;
226 
227  void computeD2Forward_() const;
228 };
229 }
230 #endif // BPP_NUMERIC_HMM_RESCALEDHMMLIKELIHOOD_H
void setBreakPoints(const std::vector< size_t > &breakPoints) override
Hidden states alphabet.
double getLogLikelihood() const override
std::vector< size_t > breakPoints_
HmmStateAlphabet & hmmStateAlphabet() override
void computeDLikelihood_() const override
std::vector< double > likelihood_
The likelihood arrays.
std::vector< double > dScales_
double getLikelihoodForASite(size_t site) const override
Get the likelihood for a site, and its derivatives.
A partial implementation of the Parametrizable interface.
HmmEmissionProbabilities & hmmEmissionProbabilities() override
Interface for computing emission probabilities in a Hidden Markov Model.
void setNamespace(const std::string &nameSpace) override
Set the namespace for the parameter names.
std::shared_ptr< const HmmEmissionProbabilities > getHmmEmissionProbabilities() const override
void getHiddenStatesPosteriorProbabilities(std::vector< std::vector< double >> &probs, bool append=false) const override
std::vector< std::vector< double > > backLikelihood_
backward likelihood
The parameter list object.
Definition: ParameterList.h:27
std::shared_ptr< HmmEmissionProbabilities > emissionProbabilities_
const HmmStateAlphabet & hmmStateAlphabet() const override
RescaledHmmLikelihood(std::shared_ptr< HmmStateAlphabet > hiddenAlphabet, std::shared_ptr< HmmTransitionMatrix > transitionMatrix, std::shared_ptr< HmmEmissionProbabilities > emissionProbabilities, const std::string &prefix)
Build a new RescaledHmmLikelihood object.
std::vector< double > d2Scales_
std::shared_ptr< HmmStateAlphabet > getHmmStateAlphabet() override
void setParametersValues(const ParameterList &parameters) override
Update the parameters from parameters.
void setParameters(const ParameterList &pl) override
Set the point where the function must be computed.
Vdouble getLikelihoodForEachSite() const override
Get the likelihood for each site.
void fireParameterChanged(const ParameterList &pl) override
Notify the class when one or several parameters have changed.
const HmmTransitionMatrix & hmmTransitionMatrix() const override
std::vector< double > Vdouble
Definition: VectorTools.h:34
std::vector< double > scales_
scales for likelihood computing
A simple implementation of hidden Markov models recursion.
RescaledHmmLikelihood(const RescaledHmmLikelihood &lik)
std::vector< std::vector< double > > d2Likelihood_
partial impmementation of Hmm Likelihoods.
std::shared_ptr< const HmmTransitionMatrix > getHmmTransitionMatrix() const override
const std::vector< size_t > & getBreakPoints() const override
RescaledHmmLikelihood & operator=(const RescaledHmmLikelihood &lik)
std::shared_ptr< HmmStateAlphabet > hiddenAlphabet_
The alphabet describing the hidden states.
AbstractHmmLikelihood & operator=(const AbstractHmmLikelihood &adhlik)
std::shared_ptr< HmmEmissionProbabilities > getHmmEmissionProbabilities() override
std::shared_ptr< HmmTransitionMatrix > transitionMatrix_
double getDLogLikelihoodForASite(size_t site) const override
RescaledHmmLikelihood * clone() const override
Create a copy of this object and send a pointer to it.
HmmTransitionMatrix & hmmTransitionMatrix() override
std::shared_ptr< HmmTransitionMatrix > getHmmTransitionMatrix() override
Describe the transition probabilities between hidden states of a Hidden Markov Model.
void computeD2Likelihood_() const override
Vdouble getHiddenStatesPosteriorProbabilitiesForASite(size_t site) const override
double getValue() const override
Get the value of the function at the current point.
std::shared_ptr< const HmmStateAlphabet > getHmmStateAlphabet() const override
const HmmEmissionProbabilities & hmmEmissionProbabilities() const override
double getD2LogLikelihoodForASite(size_t site) const override
std::vector< std::vector< double > > dLikelihood_
derivatec of forward likelihood