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