bpp-core3  3.0.0
LogsumHmmLikelihood.h
Go to the documentation of this file.
1 //
2 // File: LogsumHmmLikelihood.h
3 // Authors:
4 // Julien Dutheil
5 // Created: 2009-11-23 11:27: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_LOGSUMHMMLIKELIHOOD_H
42 #define BPP_NUMERIC_HMM_LOGSUMHMMLIKELIHOOD_H
43 
44 
45 #include "../AbstractParametrizable.h"
46 #include "../Matrix/Matrix.h"
47 #include "../NumTools.h"
48 #include "HmmLikelihood.h"
49 
50 // From the STL:
51 #include <vector>
52 #include <memory>
53 
54 namespace bpp
55 {
68  public virtual AbstractHmmLikelihood,
70 {
71 protected:
79 
86 
93  std::vector<double> logLikelihood_;
94  std::vector<double> partialLogLikelihoods_;
95  double logLik_;
96 
108  mutable std::vector< std::vector<double> > dLogLikelihood_;
109  mutable std::vector<double> partialDLogLikelihoods_;
110 
111  mutable std::vector< std::vector<double> > d2LogLikelihood_;
112  mutable std::vector<double> partialD2LogLikelihoods_;
113 
121  mutable std::vector<std::vector<double> > backLogLikelihood_;
123 
124  std::vector<size_t> breakPoints_;
125 
127 
128 public:
143  HmmStateAlphabet* hiddenAlphabet,
144  HmmTransitionMatrix* transitionMatrix,
145  HmmEmissionProbabilities* emissionProbabilities,
146  bool ownsPointers_ = true,
147  const std::string& prefix = "");
148 
152  hiddenAlphabet_(),
158  logLik_(lik.logLik_),
166  nbStates_(lik.nbStates_),
167  nbSites_(lik.nbSites_)
168  {
169  if (ownsPointers_)
170  {
171  hiddenAlphabet_ = dynamic_cast<HmmStateAlphabet*>(lik.hiddenAlphabet_->clone());
174  }
175  else
176  {
180  }
181 
182  // Now adjust pointers:
185  }
186 
188  {
190  AbstractParametrizable::operator=(lik);
191 
193 
194  if (ownsPointers_)
195  {
196  hiddenAlphabet_ = dynamic_cast<HmmStateAlphabet*>(lik.hiddenAlphabet_->clone());
199  }
200  else
201  {
205  }
206 
215  logLik_ = lik.logLik_;
217  nbStates_ = lik.nbStates_;
218  nbSites_ = lik.nbSites_;
219 
220  // Now adjust pointers:
223  return *this;
224  }
225 
227  {
228  if (ownsPointers_)
229  {
230  delete hiddenAlphabet_;
231  delete transitionMatrix_;
232  delete emissionProbabilities_;
233  }
234  }
235 
236 
237  LogsumHmmLikelihood* clone() const { return new LogsumHmmLikelihood(*this); }
238 
239 public:
242 
245 
248 
249  void setBreakPoints(const std::vector<size_t>& breakPoints)
250  {
251  breakPoints_ = breakPoints;
252  computeForward_();
254  }
255 
256  const std::vector<size_t>& getBreakPoints() const { return breakPoints_; }
257 
258  void setParameters(const ParameterList& pl)
259  {
261  }
262 
263  double getValue() const { return -logLik_; }
264 
265  double getLogLikelihood() const { return logLik_; }
266 
267  double getLikelihoodForASite(size_t site) const;
268 
269  double getDLogLikelihoodForASite(size_t site) const;
270 
271  double getD2LogLikelihoodForASite(size_t site) const;
272 
274 
275  void setNamespace(const std::string& nameSpace);
276 
277  void fireParameterChanged(const ParameterList& pl);
278 
280 
281  void getHiddenStatesPosteriorProbabilities(std::vector< std::vector<double> >& probs, bool append = false) const;
282 
283  void computeLikelihood();
284 
285 protected:
286  void computeForward_();
287 
288  void computeBackward_() const;
289 
290  void computeDLikelihood_() const
291  {
293  }
294 
295  void computeD2Likelihood_() const
296  {
298  }
299 
300  void computeDForward_() const;
301 
302  void computeD2Forward_() const;
303 };
304 }
305 #endif // BPP_NUMERIC_HMM_LOGSUMHMMLIKELIHOOD_H
AbstractHmmLikelihood & operator=(const AbstractHmmLikelihood &adhlik)
A partial implementation of the Parametrizable interface.
void setParametersValues(const ParameterList &parameters)
Update the parameters from parameters.
virtual Clonable * clone() const =0
Create a copy of this object and send a pointer to it.
Interface for computing emission probabilities in a Hidden Markov Model.
virtual HmmEmissionProbabilities * clone() const =0
Create a copy of this object and send a pointer to it.
virtual void setHmmStateAlphabet(const HmmStateAlphabet *stateAlphabet)=0
Set the new hidden state alphabet.
Hidden states alphabet.
Describe the transition probabilities between hidden states of a Hidden Markov Model.
virtual void setHmmStateAlphabet(const HmmStateAlphabet *stateAlphabet)=0
Set the new hidden state alphabet.
A simple implementation of hidden Markov models recursion.
HmmStateAlphabet * hiddenAlphabet_
The alphabet describing the hidden states.
const HmmStateAlphabet & getHmmStateAlphabet() const
const HmmTransitionMatrix & getHmmTransitionMatrix() const
double getLikelihoodForASite(size_t site) const
Get the likelihood for a site, and its derivatives.
std::vector< size_t > breakPoints_
std::vector< double > partialD2LogLikelihoods_
Vdouble getLikelihoodForEachSite() const
Get the likelihood for each site.
HmmTransitionMatrix & getHmmTransitionMatrix()
LogsumHmmLikelihood * clone() const
Create a copy of this object and send a pointer to it.
std::vector< double > partialDLogLikelihoods_
double getDLogLikelihoodForASite(size_t site) const
HmmStateAlphabet & getHmmStateAlphabet()
void setNamespace(const std::string &nameSpace)
Set the namespace for the parameter names.
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.
HmmEmissionProbabilities * emissionProbabilities_
LogsumHmmLikelihood(const LogsumHmmLikelihood &lik)
std::vector< std::vector< double > > backLogLikelihood_
backward logLikelihood
double getValue() const
Get the value of the function at the current point.
bool ownsPointers_
Owns previous objects.
std::vector< double > partialLogLikelihoods_
void setBreakPoints(const std::vector< size_t > &breakPoints)
std::vector< std::vector< double > > d2LogLikelihood_
HmmEmissionProbabilities & getHmmEmissionProbabilities()
LogsumHmmLikelihood(HmmStateAlphabet *hiddenAlphabet, HmmTransitionMatrix *transitionMatrix, HmmEmissionProbabilities *emissionProbabilities, bool ownsPointers_=true, const std::string &prefix="")
Build a new LogsumHmmLikelihood object.
const HmmEmissionProbabilities & getHmmEmissionProbabilities() const
HmmTransitionMatrix * transitionMatrix_
Vdouble getHiddenStatesPosteriorProbabilitiesForASite(size_t site) const
void setParameters(const ParameterList &pl)
Set the point where the function must be computed.
LogsumHmmLikelihood & operator=(const LogsumHmmLikelihood &lik)
std::vector< double > logLikelihood_
The likelihood array.
std::vector< std::vector< double > > dLogLikelihood_
The DLogLikelihood arrays.
const std::vector< size_t > & getBreakPoints() const
double getD2LogLikelihoodForASite(size_t site) const
The parameter list object.
Definition: ParameterList.h:65
std::vector< double > Vdouble
Definition: VectorTools.h:70