bpp-core  2.2.0
RescaledHmmLikelihood.h
Go to the documentation of this file.
1 //
2 // File: RescaledHmmLikelihood.h
3 // Created by: Julien Dutheil
4 // Created on: Fri Oct 26 11:57 2007
5 //
6 
7 /*
8  Copyright or © or Copr. CNRS, (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 _RESCALEDHMMLIKELIHOOD_H_
41 #define _RESCALEDHMMLIKELIHOOD_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 
60  public virtual AbstractHmmLikelihood,
62  {
63  private:
67  std::auto_ptr<HmmStateAlphabet> hiddenAlphabet_;
68  std::auto_ptr<HmmTransitionMatrix> transitionMatrix_;
69  std::auto_ptr<HmmEmissionProbabilities> emissionProbabilities_;
70 
83  std::vector<double> likelihood_;
84 
92  mutable std::vector<std::vector<double> > dLikelihood_;
93  mutable std::vector<std::vector<double> > d2Likelihood_;
94 
102  mutable std::vector<std::vector<double> > backLikelihood_;
104 
112  std::vector<double> scales_;
113 
114  mutable std::vector<double> dScales_;
115  mutable std::vector<double> d2Scales_;
116  double logLik_;
117 
118  std::vector<size_t> breakPoints_;
119 
121 
122  public:
132  HmmStateAlphabet* hiddenAlphabet,
133  HmmTransitionMatrix* transitionMatrix,
134  HmmEmissionProbabilities* emissionProbabilities,
135  const std::string& prefix) throw (Exception);
136 
140  hiddenAlphabet_(dynamic_cast<HmmStateAlphabet*>(lik.hiddenAlphabet_->clone())),
148  scales_(lik.scales_),
149  dScales_(lik.dScales_),
150  d2Scales_(lik.d2Scales_),
151  logLik_(lik.logLik_),
153  nbStates_(lik.nbStates_),
154  nbSites_(lik.nbSites_)
155  {
156  // Now adjust pointers:
157  transitionMatrix_->setHmmStateAlphabet(hiddenAlphabet_.get());
158  emissionProbabilities_->setHmmStateAlphabet(hiddenAlphabet_.get());
159  }
160 
162  {
164  AbstractParametrizable::operator=(lik);
165  hiddenAlphabet_ = std::auto_ptr<HmmStateAlphabet>(dynamic_cast<HmmStateAlphabet*>(lik.hiddenAlphabet_->clone()));
166  transitionMatrix_ = std::auto_ptr<HmmTransitionMatrix>(dynamic_cast<HmmTransitionMatrix*>(lik.transitionMatrix_->clone()));
167  emissionProbabilities_ = std::auto_ptr<HmmEmissionProbabilities>(dynamic_cast<HmmEmissionProbabilities*>(lik.emissionProbabilities_->clone()));
168  likelihood_ = lik.likelihood_;
173  scales_ = lik.scales_;
174  dScales_ = lik.dScales_;
175  d2Scales_ = lik.d2Scales_;
176  logLik_ = lik.logLik_;
178  nbStates_ = lik.nbStates_;
179  nbSites_ = lik.nbSites_;
180 
181  // Now adjust pointers:
182  transitionMatrix_->setHmmStateAlphabet(hiddenAlphabet_.get());
183  emissionProbabilities_->setHmmStateAlphabet(hiddenAlphabet_.get());
184  return *this;
185  }
186 
188 
189 #ifndef NO_VIRTUAL_COV
191 #else
192  Clonable*
193 #endif
194  clone() const { return new RescaledHmmLikelihood(*this); }
195 
196  public:
199 
202 
205 
206  void setBreakPoints(const std::vector<size_t>& breakPoints) {
207  breakPoints_ = breakPoints;
208  computeForward_();
210  }
211 
212  const std::vector<size_t>& getBreakPoints() const { return breakPoints_; }
213 
214  void setParameters(const ParameterList& pl) throw (Exception)
215  {
217  }
218 
219  double getValue() const throw (Exception) { return -logLik_; }
220 
221  double getLogLikelihood() const { return logLik_; }
222 
223  double getLikelihoodForASite(size_t site) const;
224 
226 
227  void fireParameterChanged(const ParameterList& pl);
228 
230 
231  void getHiddenStatesPosteriorProbabilities(std::vector< std::vector<double> >& probs, bool append = false) const throw (Exception);
232 
233 
234  protected:
235  void computeForward_();
236  void computeBackward_() const;
237 
238 
239  void computeDLikelihood_() const
240  {
242  }
243 
244  void computeD2Likelihood_() const
245  {
247  }
248 
249  void computeDForward_() const;
250 
251  void computeD2Forward_() const;
252 
253  };
254 
255 }
256 
257 #endif //_RESCALEDHMMLIKELIHOOD_H_
258 
Hidden states alphabet.
std::vector< size_t > breakPoints_
const HmmStateAlphabet & getHmmStateAlphabet() const
std::vector< double > likelihood_
The likelihood arrays.
HmmEmissionProbabilities & getHmmEmissionProbabilities()
This class allows to perform a correspondence analysis.
void fireParameterChanged(const ParameterList &pl)
Notify the class when one or several parameters have changed.
std::vector< double > dScales_
A partial implementation of the Parametrizable interface.
const std::vector< size_t > & getBreakPoints() const
void setBreakPoints(const std::vector< size_t > &breakPoints)
Vdouble getHiddenStatesPosteriorProbabilitiesForASite(size_t site) const
RescaledHmmLikelihood * clone() const
Create a copy of this object and send a pointer to it.
const HmmEmissionProbabilities & getHmmEmissionProbabilities() const
Interface for computing emission probabilities in a Hidden Markov Model.
std::auto_ptr< HmmTransitionMatrix > transitionMatrix_
std::vector< std::vector< double > > backLikelihood_
backward likelihood
void setParameters(const ParameterList &pl)
Set the point where the function must be computed.
The parameter list object.
Definition: ParameterList.h:61
RescaledHmmLikelihood(HmmStateAlphabet *hiddenAlphabet, HmmTransitionMatrix *transitionMatrix, HmmEmissionProbabilities *emissionProbabilities, const std::string &prefix)
Build a new RescaledHmmLikelihood object.
void getHiddenStatesPosteriorProbabilities(std::vector< std::vector< double > > &probs, bool append=false) const
std::vector< double > d2Scales_
std::vector< double > Vdouble
Definition: VectorTools.h:67
std::vector< double > scales_
scales for likelihood computing
Vdouble getLikelihoodForEachSite() const
Get the likelihood for each site.
A simple implementation of hidden Markov models recursion.
RescaledHmmLikelihood(const RescaledHmmLikelihood &lik)
std::vector< std::vector< double > > d2Likelihood_
void setParametersValues(const ParameterList &parameters)
Update the parameters from parameters.
RescaledHmmLikelihood & operator=(const RescaledHmmLikelihood &lik)
Exception base class.
Definition: Exceptions.h:57
The Clonable interface (allow an object to be cloned).
Definition: Clonable.h:99
HmmTransitionMatrix & getHmmTransitionMatrix()
AbstractHmmLikelihood & operator=(const AbstractHmmLikelihood &adhlik)
HmmStateAlphabet & getHmmStateAlphabet()
double getValue() const
Get the value of the function at the current point.
std::auto_ptr< HmmEmissionProbabilities > emissionProbabilities_
double getLikelihoodForASite(size_t site) const
Get the likelihood for a site.
Describe the transition probabilities between hidden states of a Hidden Markov Model.
const HmmTransitionMatrix & getHmmTransitionMatrix() const
std::auto_ptr< HmmStateAlphabet > hiddenAlphabet_
The alphabet describing the hidden states.
std::vector< std::vector< double > > dLikelihood_
derivatec of forward likelihood