bpp-core3  3.0.0
LowMemoryRescaledHmmLikelihood.h
Go to the documentation of this file.
1 //
2 // File: LowMemoryRescaledHmmLikelihood.h
3 // Authors:
4 // Julien Dutheil
5 // Created: 2009-12-16 10:47:00
6 //
7 
8 /*
9  Copyright or © or Copr. Bio++ Development Team, (November 16, 2004)
10 
11  This software is a computer program whose purpose is to provide classes
12  for phylogenetic data analysis.
13 
14  This software is governed by the CeCILL license under French law and
15  abiding by the rules of distribution of free software. You can use,
16  modify and/ or redistribute the software under the terms of the CeCILL
17  license as circulated by CEA, CNRS and INRIA at the following URL
18  "http://www.cecill.info".
19 
20  As a counterpart to the access to the source code and rights to copy,
21  modify and redistribute granted by the license, users are provided only
22  with a limited warranty and the software's author, the holder of the
23  economic rights, and the successive licensors have only limited
24  liability.
25 
26  In this respect, the user's attention is drawn to the risks associated
27  with loading, using, modifying and/or developing or reproducing the
28  software by the user in light of its specific status of free software,
29  that may mean that it is complicated to manipulate, and that also
30  therefore means that it is reserved for developers and experienced
31  professionals having in-depth computer knowledge. Users are therefore
32  encouraged to load and test the software's suitability as regards their
33  requirements in conditions enabling the security of their systems and/or
34  data to be ensured and, more generally, to use and operate it in the
35  same conditions as regards security.
36 
37  The fact that you are presently reading this means that you have had
38  knowledge of the CeCILL license and that you accept its terms.
39 */
40 
41 #ifndef BPP_NUMERIC_HMM_LOWMEMORYRESCALEDHMMLIKELIHOOD_H
42 #define BPP_NUMERIC_HMM_LOWMEMORYRESCALEDHMMLIKELIHOOD_H
43 
44 
45 #include "../AbstractParametrizable.h"
46 #include "../Matrix/Matrix.h"
47 #include "HmmLikelihood.h"
48 
49 // From the STL:
50 #include <vector>
51 #include <memory>
52 
53 namespace bpp
54 {
69  public AbstractHmmLikelihood,
71 {
72 private:
76  std::unique_ptr<HmmStateAlphabet> hiddenAlphabet_;
77  std::unique_ptr<HmmTransitionMatrix> transitionMatrix_;
78  std::unique_ptr<HmmEmissionProbabilities> emissionProbabilities_;
79 
85  std::vector<double> likelihood1_;
86  std::vector<double> likelihood2_;
87  double logLik_;
88  size_t maxSize_;
89 
90  std::vector<size_t> breakPoints_;
91 
93 
94 public:
111  HmmStateAlphabet* hiddenAlphabet,
112  HmmTransitionMatrix* transitionMatrix,
113  HmmEmissionProbabilities* emissionProbabilities,
114  const std::string& prefix,
115  size_t maxSize = 1000000);
116 
120  hiddenAlphabet_(dynamic_cast<HmmStateAlphabet*>(lik.hiddenAlphabet_->clone())),
125  logLik_(lik.logLik_),
126  maxSize_(lik.maxSize_),
128  nbStates_(lik.nbStates_),
129  nbSites_(lik.nbSites_)
130  {
131  // Now adjust pointers:
132  transitionMatrix_->setHmmStateAlphabet(hiddenAlphabet_.get());
133  emissionProbabilities_->setHmmStateAlphabet(hiddenAlphabet_.get());
134  }
135 
137  {
139  AbstractParametrizable::operator=(lik);
140  hiddenAlphabet_ = std::unique_ptr<HmmStateAlphabet>(dynamic_cast<HmmStateAlphabet*>(lik.hiddenAlphabet_->clone()));
141  transitionMatrix_ = std::unique_ptr<HmmTransitionMatrix>(dynamic_cast<HmmTransitionMatrix*>(lik.transitionMatrix_->clone()));
142  emissionProbabilities_ = std::unique_ptr<HmmEmissionProbabilities>(dynamic_cast<HmmEmissionProbabilities*>(lik.emissionProbabilities_->clone()));
145  logLik_ = lik.logLik_;
146  maxSize_ = lik.maxSize_;
148  nbStates_ = lik.nbStates_;
149  nbSites_ = lik.nbSites_;
150 
151  // Now adjust pointers:
152  transitionMatrix_->setHmmStateAlphabet(hiddenAlphabet_.get());
153  emissionProbabilities_->setHmmStateAlphabet(hiddenAlphabet_.get());
154  return *this;
155  }
156 
158 
160 
161 public:
164 
167 
170 
171  void setBreakPoints(const std::vector<size_t>& breakPoints)
172  {
173  breakPoints_ = breakPoints;
174  computeForward_();
175  }
176 
177  const std::vector<size_t>& getBreakPoints() const { return breakPoints_; }
178 
179  void setParameters(const ParameterList& pl)
180  {
182  }
183 
184  double getValue() const { return -logLik_; }
185 
186  double getLogLikelihood() const { return logLik_; }
187 
188  void setNamespace(const std::string& nameSpace);
189 
190  void fireParameterChanged(const ParameterList& pl);
191 
192  double getLikelihoodForASite(size_t site) const
193  {
194  throw (NotImplementedException("LowMemoryRescaledHmmLikelihood::getLikelihoodForASite. This class can't compute posterior probabilities, use RescaledHmmLikelihood instead."));
195  }
196 
197 
199  {
200  throw (NotImplementedException("LowMemoryRescaledHmmLikelihood::getLikelihoodForEachSite. This class can't compute posterior probabilities, use RescaledHmmLikelihood instead."));
201  }
202 
204  {
205  throw (NotImplementedException("LowMemoryRescaledHmmLikelihood::getHiddenStatesPosteriorProbabilitiesForASite. This class can't compute posterior probabilities, use RescaledHmmLikelihood instead."));
206  }
207 
208  void getHiddenStatesPosteriorProbabilities(std::vector< std::vector<double> >& probs, bool append = false) const
209  {
210  throw (NotImplementedException("LowMemoryRescaledHmmLikelihood::getHiddenStatesPosteriorProbabilities. This class can't compute posterior probabilities, use RescaledHmmLikelihood instead."));
211  }
212 
213 protected:
214  void computeForward_();
215 
216  void computeDLikelihood_() const
217  {
218  throw (NotImplementedException("LowMemoryRescaledHmmLikelihood::computeDLikelihood_. Use RescaledHmmLikelihood instead."));
219  }
220 
221  void computeD2Likelihood_() const
222  {
223  throw (NotImplementedException("LowMemoryRescaledHmmLikelihood::computeD2Likelihood_. Use RescaledHmmLikelihood instead."));
224  }
225 
226  double getDLogLikelihoodForASite(size_t site) const
227  {
228  throw (NotImplementedException("LowMemoryRescaledHmmLikelihood::getDLogLikelihoodForASite. Use RescaledHmmLikelihood instead."));
229  return 0;
230  }
231 
232  double getD2LogLikelihoodForASite(size_t site) const
233  {
234  throw (NotImplementedException("LowMemoryRescaledHmmLikelihood::getD2LogLikelihoodForASite. Use RescaledHmmLikelihood instead."));
235  return 0;
236  }
237 };
238 }
239 #endif // BPP_NUMERIC_HMM_LOWMEMORYRESCALEDHMMLIKELIHOOD_H
AbstractHmmLikelihood & operator=(const AbstractHmmLikelihood &adhlik)
A partial implementation of the Parametrizable interface.
void setParametersValues(const ParameterList &parameters)
Update the parameters from parameters.
Interface for computing emission probabilities in a Hidden Markov Model.
Hidden states alphabet.
Describe the transition probabilities between hidden states of a Hidden Markov Model.
A modified implementation of the RescaledHmmLikelihood implementation, with lower memory usage.
const HmmTransitionMatrix & getHmmTransitionMatrix() const
void setNamespace(const std::string &nameSpace)
Set the namespace for the parameter names.
LowMemoryRescaledHmmLikelihood * clone() const
Create a copy of this object and send a pointer to it.
std::unique_ptr< HmmStateAlphabet > hiddenAlphabet_
The alphabet describing the hidden states.
const HmmEmissionProbabilities & getHmmEmissionProbabilities() const
Vdouble getHiddenStatesPosteriorProbabilitiesForASite(size_t site) const
const std::vector< size_t > & getBreakPoints() const
Vdouble getLikelihoodForEachSite() const
Get the likelihood for each site.
double getLikelihoodForASite(size_t site) const
Get the likelihood for a site, and its derivatives.
std::unique_ptr< HmmTransitionMatrix > transitionMatrix_
LowMemoryRescaledHmmLikelihood & operator=(const LowMemoryRescaledHmmLikelihood &lik)
void setBreakPoints(const std::vector< size_t > &breakPoints)
HmmEmissionProbabilities & getHmmEmissionProbabilities()
const HmmStateAlphabet & getHmmStateAlphabet() const
void getHiddenStatesPosteriorProbabilities(std::vector< std::vector< double > > &probs, bool append=false) const
void fireParameterChanged(const ParameterList &pl)
Notify the class when one or several parameters have changed.
void setParameters(const ParameterList &pl)
Set the point where the function must be computed.
double getValue() const
Get the value of the function at the current point.
std::vector< double > likelihood1_
The likelihood array.
LowMemoryRescaledHmmLikelihood(const LowMemoryRescaledHmmLikelihood &lik)
LowMemoryRescaledHmmLikelihood(HmmStateAlphabet *hiddenAlphabet, HmmTransitionMatrix *transitionMatrix, HmmEmissionProbabilities *emissionProbabilities, const std::string &prefix, size_t maxSize=1000000)
Build a new LowMemoryRescaledHmmLikelihood object.
std::unique_ptr< HmmEmissionProbabilities > emissionProbabilities_
This exception is sent when a given method is not implemented.
Definition: Exceptions.h:230
The parameter list object.
Definition: ParameterList.h:65
std::vector< double > Vdouble
Definition: VectorTools.h:70