bpp-core3  3.0.0
HmmLikelihood.h
Go to the documentation of this file.
1 //
2 // File: HmmLikelihood.h
3 // Authors:
4 // Julien Dutheil
5 // Created: 2007-10-26 11:20: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_HMMLIKELIHOOD_H
42 #define BPP_NUMERIC_HMM_HMMLIKELIHOOD_H
43 
44 
45 
46 // From NumCalc:
47 #include "../Function/Functions.h"
48 #include "../VectorTools.h"
49 
50 #include "HmmStateAlphabet.h"
51 #include "HmmTransitionMatrix.h"
53 
54 namespace bpp
55 {
79  public virtual DerivableSecondOrder
80 {
81 public:
82  virtual HmmLikelihood* clone() const = 0;
83 
84  virtual const HmmStateAlphabet& getHmmStateAlphabet() const = 0;
86 
87  virtual const HmmTransitionMatrix& getHmmTransitionMatrix() const = 0;
89 
92 
93  virtual void getHiddenStatesPosteriorProbabilities(std::vector< std::vector<double> >& probs, bool append) const = 0;
94 
96 
97  virtual double getLogLikelihood() const = 0;
98 
99  virtual double getDLogLikelihood() const = 0;
100 
101  virtual double getD2LogLikelihood() const = 0;
102 
110  virtual double getLikelihoodForASite(size_t site) const = 0;
111 
112  virtual double getDLogLikelihoodForASite(size_t site) const = 0;
113 
114  virtual double getD2LogLikelihoodForASite(size_t site) const = 0;
120  virtual Vdouble getLikelihoodForEachSite() const = 0;
121 
122  virtual const std::vector<size_t>& getBreakPoints() const = 0;
123 
124  virtual void setBreakPoints(const std::vector<size_t>& breakPoints) = 0;
125 
126 protected:
127  virtual void computeDLikelihood_() const = 0;
128 
129  virtual void computeD2Likelihood_() const = 0;
130 };
131 
132 
133 /*
134  * @brief partial impmementation of Hmm Likelihoods.
135  *
136  */
137 
139  public virtual HmmLikelihood
140 {
141 protected:
142  mutable double dLogLik_;
143  mutable std::string dVariable_;
144 
145  mutable double d2LogLik_;
146  mutable std::string d2Variable_;
147 
148 public:
150 
152 
154 
155  /* @{
156  *
157  * @brief From FirstOrder:
158  *
159  */
161 
162  bool enableFirstOrderDerivatives() const { return true;}
163 
164  double getFirstOrderDerivative(const std::string& variable) const;
165 
166  double getDLogLikelihood() const
167  {
168  return dLogLik_;
169  }
170 
171  /*
172  * @}
173  *
174  */
175 
176  /* @{
177  *
178  * @brief From SecondOrder:
179  *
180  */
182 
183  bool enableSecondOrderDerivatives() const {return true;}
184 
185  double getSecondOrderDerivative(const std::string& variable) const;
186 
187  double getD2LogLikelihood() const
188  {
189  return d2LogLik_;
190  }
191 
192  double getSecondOrderDerivative(const std::string& variable1, const std::string& variable2) const
193  {
194  throw (NotImplementedException("AbstractHmmLikelihood::getSecondOrderDerivative is not defined for 2 variables."));
195  }
196 
197  /*
198  * @}
199  *
200  */
201 };
202 } // end of namespace bpp.
203 #endif // BPP_NUMERIC_HMM_HMMLIKELIHOOD_H
double getSecondOrderDerivative(const std::string &variable) const
Get the second order derivative of the function at the current point.
double getD2LogLikelihood() const
AbstractHmmLikelihood & operator=(const AbstractHmmLikelihood &adhlik)
bool enableSecondOrderDerivatives() const
Tell if derivatives must be computed.
void enableSecondOrderDerivatives(bool yn)
Tell if derivatives must be computed.
void enableFirstOrderDerivatives(bool yn)
Tell if derivatives must be computed.
double getDLogLikelihood() const
bool enableFirstOrderDerivatives() const
Tell if derivatives must be computed.
double getFirstOrderDerivative(const std::string &variable) const
Get the derivative of the function at the current point.
double getSecondOrderDerivative(const std::string &variable1, const std::string &variable2) const
Get the value of the cross derivative of the function according to a given set of parameters.
This is the abstract class for second order derivable functions.
Definition: Functions.h:188
Interface for computing emission probabilities in a Hidden Markov Model.
Basal interface for Hidden Markov Models likelihood computation.
Definition: HmmLikelihood.h:80
virtual void computeD2Likelihood_() const =0
virtual const std::vector< size_t > & getBreakPoints() const =0
virtual double getLikelihoodForASite(size_t site) const =0
Get the likelihood for a site, and its derivatives.
virtual double getD2LogLikelihoodForASite(size_t site) const =0
virtual void computeDLikelihood_() const =0
virtual double getDLogLikelihoodForASite(size_t site) const =0
virtual double getLogLikelihood() const =0
virtual HmmEmissionProbabilities & getHmmEmissionProbabilities()=0
virtual double getD2LogLikelihood() const =0
virtual double getDLogLikelihood() const =0
virtual Vdouble getHiddenStatesPosteriorProbabilitiesForASite(size_t site) const =0
virtual HmmStateAlphabet & getHmmStateAlphabet()=0
virtual const HmmStateAlphabet & getHmmStateAlphabet() const =0
virtual const HmmEmissionProbabilities & getHmmEmissionProbabilities() const =0
virtual void getHiddenStatesPosteriorProbabilities(std::vector< std::vector< double > > &probs, bool append) const =0
virtual Vdouble getLikelihoodForEachSite() const =0
Get the likelihood for each site.
virtual HmmTransitionMatrix & getHmmTransitionMatrix()=0
virtual const HmmTransitionMatrix & getHmmTransitionMatrix() const =0
virtual void setBreakPoints(const std::vector< size_t > &breakPoints)=0
virtual HmmLikelihood * clone() const =0
Create a copy of this object and send a pointer to it.
Hidden states alphabet.
Describe the transition probabilities between hidden states of a Hidden Markov Model.
This exception is sent when a given method is not implemented.
Definition: Exceptions.h:230
std::vector< double > Vdouble
Definition: VectorTools.h:70