bpp-core3  3.0.0
LowMemoryRescaledHmmLikelihood.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_LOWMEMORYRESCALEDHMMLIKELIHOOD_H
6 #define BPP_NUMERIC_HMM_LOWMEMORYRESCALEDHMMLIKELIHOOD_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 {
33  public AbstractHmmLikelihood,
35 {
36 private:
40  std::shared_ptr<HmmStateAlphabet> hiddenAlphabet_;
41  std::shared_ptr<HmmTransitionMatrix> transitionMatrix_;
42  std::shared_ptr<HmmEmissionProbabilities> emissionProbabilities_;
43 
49  std::vector<double> likelihood1_;
50  std::vector<double> likelihood2_;
51  double logLik_;
52  size_t maxSize_;
53 
54  std::vector<size_t> breakPoints_;
55 
57 
58 public:
75  std::shared_ptr<HmmStateAlphabet> hiddenAlphabet,
76  std::shared_ptr<HmmTransitionMatrix> transitionMatrix,
77  std::shared_ptr<HmmEmissionProbabilities> emissionProbabilities,
78  const std::string& prefix,
79  size_t maxSize = 1000000);
80 
84  hiddenAlphabet_(dynamic_cast<HmmStateAlphabet*>(lik.hiddenAlphabet_->clone())),
85  transitionMatrix_(dynamic_cast<HmmTransitionMatrix*>(lik.transitionMatrix_->clone())),
86  emissionProbabilities_(dynamic_cast<HmmEmissionProbabilities*>(lik.emissionProbabilities_->clone())),
87  likelihood1_(lik.likelihood1_),
88  likelihood2_(lik.likelihood2_),
89  logLik_(lik.logLik_),
90  maxSize_(lik.maxSize_),
91  breakPoints_(lik.breakPoints_),
92  nbStates_(lik.nbStates_),
93  nbSites_(lik.nbSites_)
94  {
95  // Now adjust pointers:
96  transitionMatrix_->setHmmStateAlphabet(hiddenAlphabet_);
97  emissionProbabilities_->setHmmStateAlphabet(hiddenAlphabet_);
98  }
99 
101  {
103  AbstractParametrizable::operator=(lik);
104  hiddenAlphabet_ = std::unique_ptr<HmmStateAlphabet>(dynamic_cast<HmmStateAlphabet*>(lik.hiddenAlphabet_->clone()));
105  transitionMatrix_ = std::unique_ptr<HmmTransitionMatrix>(dynamic_cast<HmmTransitionMatrix*>(lik.transitionMatrix_->clone()));
106  emissionProbabilities_ = std::unique_ptr<HmmEmissionProbabilities>(dynamic_cast<HmmEmissionProbabilities*>(lik.emissionProbabilities_->clone()));
107  likelihood1_ = lik.likelihood1_;
108  likelihood2_ = lik.likelihood2_;
109  logLik_ = lik.logLik_;
110  maxSize_ = lik.maxSize_;
111  breakPoints_ = lik.breakPoints_;
112  nbStates_ = lik.nbStates_;
113  nbSites_ = lik.nbSites_;
114 
115  // Now adjust pointers:
116  transitionMatrix_->setHmmStateAlphabet(hiddenAlphabet_);
117  emissionProbabilities_->setHmmStateAlphabet(hiddenAlphabet_);
118  return *this;
119  }
120 
122 
123  LowMemoryRescaledHmmLikelihood* clone() const override { return new LowMemoryRescaledHmmLikelihood(*this); }
124 
125 public:
126  const HmmStateAlphabet& hmmStateAlphabet() const override { return *hiddenAlphabet_; }
127  std::shared_ptr<const HmmStateAlphabet> getHmmStateAlphabet() const override { return hiddenAlphabet_; }
128 
130  std::shared_ptr<HmmStateAlphabet> getHmmStateAlphabet() override { return hiddenAlphabet_; }
131 
132  const HmmTransitionMatrix& hmmTransitionMatrix() const override { return *transitionMatrix_; }
133  std::shared_ptr<const HmmTransitionMatrix> getHmmTransitionMatrix() const override { return transitionMatrix_; }
134 
136  std::shared_ptr<HmmTransitionMatrix> getHmmTransitionMatrix() override { return transitionMatrix_; }
137 
139  std::shared_ptr<const HmmEmissionProbabilities> getHmmEmissionProbabilities() const override { return emissionProbabilities_; }
140 
142  std::shared_ptr<HmmEmissionProbabilities> getHmmEmissionProbabilities() override { return emissionProbabilities_; }
143 
144  void setBreakPoints(const std::vector<size_t>& breakPoints) override
145  {
146  breakPoints_ = breakPoints;
147  computeForward_();
148  }
149 
150  const std::vector<size_t>& getBreakPoints() const override { return breakPoints_; }
151 
152  void setParameters(const ParameterList& pl) override
153  {
155  }
156 
157  double getValue() const override { return -logLik_; }
158 
159  double getLogLikelihood() const override { return logLik_; }
160 
161  void setNamespace(const std::string& nameSpace) override;
162 
163  void fireParameterChanged(const ParameterList& pl) override;
164 
165  double getLikelihoodForASite(size_t site) const override
166  {
167  throw (NotImplementedException("LowMemoryRescaledHmmLikelihood::getLikelihoodForASite. This class can't compute posterior probabilities, use RescaledHmmLikelihood instead."));
168  }
169 
170 
172  {
173  throw (NotImplementedException("LowMemoryRescaledHmmLikelihood::getLikelihoodForEachSite. This class can't compute posterior probabilities, use RescaledHmmLikelihood instead."));
174  }
175 
177  {
178  throw (NotImplementedException("LowMemoryRescaledHmmLikelihood::getHiddenStatesPosteriorProbabilitiesForASite. This class can't compute posterior probabilities, use RescaledHmmLikelihood instead."));
179  }
180 
181  void getHiddenStatesPosteriorProbabilities(std::vector< std::vector<double>>& probs, bool append = false) const override
182  {
183  throw (NotImplementedException("LowMemoryRescaledHmmLikelihood::getHiddenStatesPosteriorProbabilities. This class can't compute posterior probabilities, use RescaledHmmLikelihood instead."));
184  }
185 
186 protected:
187  void computeForward_();
188 
189  void computeDLikelihood_() const override
190  {
191  throw (NotImplementedException("LowMemoryRescaledHmmLikelihood::computeDLikelihood_. Use RescaledHmmLikelihood instead."));
192  }
193 
194  void computeD2Likelihood_() const override
195  {
196  throw (NotImplementedException("LowMemoryRescaledHmmLikelihood::computeD2Likelihood_. Use RescaledHmmLikelihood instead."));
197  }
198 
199  double getDLogLikelihoodForASite(size_t site) const override
200  {
201  throw (NotImplementedException("LowMemoryRescaledHmmLikelihood::getDLogLikelihoodForASite. Use RescaledHmmLikelihood instead."));
202  return 0;
203  }
204 
205  double getD2LogLikelihoodForASite(size_t site) const override
206  {
207  throw (NotImplementedException("LowMemoryRescaledHmmLikelihood::getD2LogLikelihoodForASite. Use RescaledHmmLikelihood instead."));
208  return 0;
209  }
210 };
211 }
212 #endif // BPP_NUMERIC_HMM_LOWMEMORYRESCALEDHMMLIKELIHOOD_H
Hidden states alphabet.
std::shared_ptr< HmmTransitionMatrix > getHmmTransitionMatrix() override
void setParameters(const ParameterList &pl) override
Set the point where the function must be computed.
std::shared_ptr< HmmTransitionMatrix > transitionMatrix_
std::shared_ptr< HmmEmissionProbabilities > emissionProbabilities_
std::shared_ptr< const HmmTransitionMatrix > getHmmTransitionMatrix() const override
void setNamespace(const std::string &nameSpace) override
Set the namespace for the parameter names.
A partial implementation of the Parametrizable interface.
const std::vector< size_t > & getBreakPoints() const override
HmmEmissionProbabilities & hmmEmissionProbabilities() override
A modified implementation of the RescaledHmmLikelihood implementation, with lower memory usage...
Interface for computing emission probabilities in a Hidden Markov Model.
HmmTransitionMatrix & hmmTransitionMatrix() override
void getHiddenStatesPosteriorProbabilities(std::vector< std::vector< double >> &probs, bool append=false) const override
double getLikelihoodForASite(size_t site) const override
Get the likelihood for a site, and its derivatives.
void setBreakPoints(const std::vector< size_t > &breakPoints) override
The parameter list object.
Definition: ParameterList.h:27
std::shared_ptr< HmmStateAlphabet > hiddenAlphabet_
The alphabet describing the hidden states.
Vdouble getLikelihoodForEachSite() const override
Get the likelihood for each site.
void setParametersValues(const ParameterList &parameters) override
Update the parameters from parameters.
std::shared_ptr< const HmmEmissionProbabilities > getHmmEmissionProbabilities() const override
std::shared_ptr< HmmStateAlphabet > getHmmStateAlphabet() override
This exception is sent when a given method is not implemented.
Definition: Exceptions.h:191
double getDLogLikelihoodForASite(size_t site) const override
const HmmEmissionProbabilities & hmmEmissionProbabilities() const override
double getValue() const override
Get the value of the function at the current point.
double getD2LogLikelihoodForASite(size_t site) const override
void fireParameterChanged(const ParameterList &pl) override
Notify the class when one or several parameters have changed.
std::vector< double > Vdouble
Definition: VectorTools.h:34
partial impmementation of Hmm Likelihoods.
AbstractHmmLikelihood & operator=(const AbstractHmmLikelihood &adhlik)
std::shared_ptr< HmmEmissionProbabilities > getHmmEmissionProbabilities() override
std::vector< double > likelihood1_
The likelihood array.
std::shared_ptr< const HmmStateAlphabet > getHmmStateAlphabet() const override
LowMemoryRescaledHmmLikelihood * clone() const override
Create a copy of this object and send a pointer to it.
LowMemoryRescaledHmmLikelihood & operator=(const LowMemoryRescaledHmmLikelihood &lik)
const HmmStateAlphabet & hmmStateAlphabet() const override
LowMemoryRescaledHmmLikelihood(std::shared_ptr< HmmStateAlphabet > hiddenAlphabet, std::shared_ptr< HmmTransitionMatrix > transitionMatrix, std::shared_ptr< HmmEmissionProbabilities > emissionProbabilities, const std::string &prefix, size_t maxSize=1000000)
Build a new LowMemoryRescaledHmmLikelihood object.
const HmmTransitionMatrix & hmmTransitionMatrix() const override
Vdouble getHiddenStatesPosteriorProbabilitiesForASite(size_t site) const override
Describe the transition probabilities between hidden states of a Hidden Markov Model.
LowMemoryRescaledHmmLikelihood(const LowMemoryRescaledHmmLikelihood &lik)