bpp-core  2.2.0
LogsumHmmLikelihood.h
Go to the documentation of this file.
1 //
2 // File: LogsumHmmLikelihood.h
3 // Created by: Julien Dutheil
4 // Created on: Mon Nov 23 11:27 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 _LOGSUMHMMLIKELIHOOD_H_
41 #define _LOGSUMHMMLIKELIHOOD_H_
42 
43 #include "HmmLikelihood.h"
44 #include "../AbstractParametrizable.h"
45 #include "../NumTools.h"
46 #include "../Matrix/Matrix.h"
47 
48 //From the STL:
49 #include <vector>
50 #include <memory>
51 
52 namespace bpp {
53 
66  public virtual AbstractHmmLikelihood,
68  {
69  protected:
73  std::auto_ptr<HmmStateAlphabet> hiddenAlphabet_;
74  std::auto_ptr<HmmTransitionMatrix> transitionMatrix_;
75  std::auto_ptr<HmmEmissionProbabilities> emissionProbabilities_;
76 
83  std::vector<double> logLikelihood_;
84  std::vector<double> partialLogLikelihoods_;
85  double logLik_;
86 
98  mutable std::vector< std::vector<double> > dLogLikelihood_;
99  mutable std::vector<double> partialDLogLikelihoods_;
100 
101  mutable std::vector< std::vector<double> > d2LogLikelihood_;
102  mutable std::vector<double> partialD2LogLikelihoods_;
103 
111  mutable std::vector<std::vector<double> > backLogLikelihood_;
113 
114  std::vector<size_t> breakPoints_;
115 
117 
118  public:
128  HmmStateAlphabet* hiddenAlphabet,
129  HmmTransitionMatrix* transitionMatrix,
130  HmmEmissionProbabilities* emissionProbabilities,
131  const std::string& prefix) throw (Exception);
132 
136  hiddenAlphabet_(dynamic_cast<HmmStateAlphabet*>(lik.hiddenAlphabet_->clone())),
141  logLik_(lik.logLik_),
149  nbStates_(lik.nbStates_),
150  nbSites_(lik.nbSites_)
151  {
152  // Now adjust pointers:
153  transitionMatrix_->setHmmStateAlphabet(hiddenAlphabet_.get());
154  emissionProbabilities_->setHmmStateAlphabet(hiddenAlphabet_.get());
155  }
156 
158  {
160  AbstractParametrizable::operator =(lik);
161  hiddenAlphabet_ = std::auto_ptr<HmmStateAlphabet>(dynamic_cast<HmmStateAlphabet*>(lik.hiddenAlphabet_->clone()));
162  transitionMatrix_ = std::auto_ptr<HmmTransitionMatrix>(dynamic_cast<HmmTransitionMatrix*>(lik.transitionMatrix_->clone()));
163  emissionProbabilities_ = std::auto_ptr<HmmEmissionProbabilities>(dynamic_cast<HmmEmissionProbabilities*>(lik.emissionProbabilities_->clone()));
172  logLik_ = lik.logLik_;
174  nbStates_ = lik.nbStates_;
175  nbSites_ = lik.nbSites_;
176 
177  // Now adjust pointers:
178  transitionMatrix_->setHmmStateAlphabet(hiddenAlphabet_.get());
179  emissionProbabilities_->setHmmStateAlphabet(hiddenAlphabet_.get());
180  return *this;
181  }
182 
183  virtual ~LogsumHmmLikelihood() {}
184 
185 #ifndef NO_VIRTUAL_COV
187 #else
188  Clonable*
189 #endif
190  clone() const { return new LogsumHmmLikelihood(*this); }
191 
192  public:
195 
198 
201 
202  void setBreakPoints(const std::vector<size_t>& breakPoints) {
203  breakPoints_ = breakPoints;
204  computeForward_();
206  }
207 
208  const std::vector<size_t>& getBreakPoints() const { return breakPoints_; }
209 
210  void setParameters(const ParameterList& pl) throw (Exception)
211  {
213  }
214 
215  double getValue() const throw (Exception) { return -logLik_; }
216 
217  double getLogLikelihood() const { return logLik_; }
218 
219  double getLikelihoodForASite(size_t site) const;
220 
222 
223  virtual void fireParameterChanged(const ParameterList& pl);
224 
226 
227  void getHiddenStatesPosteriorProbabilities(std::vector< std::vector<double> >& probs, bool append = false) const throw (Exception);
228 
229  protected:
230  void computeForward_();
231  void computeBackward_() const;
232 
233  void computeDLikelihood_() const
234  {
236  }
237 
238  void computeD2Likelihood_() const
239  {
241  }
242 
243  void computeDForward_() const;
244 
245  void computeD2Forward_() const;
246 
247 
248  };
249 
250 }
251 
252 #endif //_LOGSUMHMMLIKELIHOOD_H_
253 
Hidden states alphabet.
std::vector< double > partialDLogLikelihoods_
std::vector< double > logLikelihood_
The likelihood array.
const std::vector< size_t > & getBreakPoints() const
This class allows to perform a correspondence analysis.
A simple implementation of hidden Markov models recursion.
double getLikelihoodForASite(size_t site) const
Get the likelihood for a site.
A partial implementation of the Parametrizable interface.
LogsumHmmLikelihood * clone() const
Create a copy of this object and send a pointer to it.
Interface for computing emission probabilities in a Hidden Markov Model.
HmmStateAlphabet & getHmmStateAlphabet()
std::vector< std::vector< double > > d2LogLikelihood_
Vdouble getLikelihoodForEachSite() const
Get the likelihood for each site.
virtual void fireParameterChanged(const ParameterList &pl)
Notify the class when one or several parameters have changed.
void setBreakPoints(const std::vector< size_t > &breakPoints)
The parameter list object.
Definition: ParameterList.h:61
HmmEmissionProbabilities & getHmmEmissionProbabilities()
LogsumHmmLikelihood(HmmStateAlphabet *hiddenAlphabet, HmmTransitionMatrix *transitionMatrix, HmmEmissionProbabilities *emissionProbabilities, const std::string &prefix)
Build a new LogsumHmmLikelihood object.
std::auto_ptr< HmmTransitionMatrix > transitionMatrix_
std::vector< std::vector< double > > backLogLikelihood_
backward logLikelihood
std::vector< double > partialLogLikelihoods_
const HmmEmissionProbabilities & getHmmEmissionProbabilities() const
std::vector< double > Vdouble
Definition: VectorTools.h:67
std::auto_ptr< HmmStateAlphabet > hiddenAlphabet_
The alphabet describing the hidden states.
std::vector< double > partialD2LogLikelihoods_
const HmmStateAlphabet & getHmmStateAlphabet() const
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 getHiddenStatesPosteriorProbabilities(std::vector< std::vector< double > > &probs, bool append=false) const
AbstractHmmLikelihood & operator=(const AbstractHmmLikelihood &adhlik)
const HmmTransitionMatrix & getHmmTransitionMatrix() const
void setParameters(const ParameterList &pl)
Set the point where the function must be computed.
Vdouble getHiddenStatesPosteriorProbabilitiesForASite(size_t site) const
std::vector< std::vector< double > > dLogLikelihood_
The DLogLikelihood arrays.
LogsumHmmLikelihood & operator=(const LogsumHmmLikelihood &lik)
double getValue() const
Get the value of the function at the current point.
std::vector< size_t > breakPoints_
HmmTransitionMatrix & getHmmTransitionMatrix()
Describe the transition probabilities between hidden states of a Hidden Markov Model.
std::auto_ptr< HmmEmissionProbabilities > emissionProbabilities_
LogsumHmmLikelihood(const LogsumHmmLikelihood &lik)