bpp-core3  3.0.0
RescaledHmmLikelihood.h
Go to the documentation of this file.
1 //
2 // File: RescaledHmmLikelihood.h
3 // Authors:
4 // Julien Dutheil
5 // Created: 2007-10-26 11:57:00
6 //
7 
8 /*
9  Copyright or © or Copr. CNRS, (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_RESCALEDHMMLIKELIHOOD_H
42 #define BPP_NUMERIC_HMM_RESCALEDHMMLIKELIHOOD_H
43 
44 
45 #include "../AbstractParametrizable.h"
46 #include "../Matrix/Matrix.h"
47 #include "HmmLikelihood.h"
48 
49 // From the STL:
50 #include <vector>
51 #include <memory>
52 
53 namespace bpp
54 {
62  public virtual AbstractHmmLikelihood,
64 {
65 private:
69  std::unique_ptr<HmmStateAlphabet> hiddenAlphabet_;
70  std::unique_ptr<HmmTransitionMatrix> transitionMatrix_;
71  std::unique_ptr<HmmEmissionProbabilities> emissionProbabilities_;
72 
85  std::vector<double> likelihood_;
86 
94  mutable std::vector<std::vector<double> > dLikelihood_;
95  mutable std::vector<std::vector<double> > d2Likelihood_;
96 
104  mutable std::vector<std::vector<double> > backLikelihood_;
106 
114  std::vector<double> scales_;
115 
116  mutable std::vector<double> dScales_;
117  mutable std::vector<double> d2Scales_;
118  double logLik_;
119 
120  std::vector<size_t> breakPoints_;
121 
123 
124 public:
134  HmmStateAlphabet* hiddenAlphabet,
135  HmmTransitionMatrix* transitionMatrix,
136  HmmEmissionProbabilities* emissionProbabilities,
137  const std::string& prefix);
138 
142  hiddenAlphabet_(dynamic_cast<HmmStateAlphabet*>(lik.hiddenAlphabet_->clone())),
150  scales_(lik.scales_),
151  dScales_(lik.dScales_),
152  d2Scales_(lik.d2Scales_),
153  logLik_(lik.logLik_),
155  nbStates_(lik.nbStates_),
156  nbSites_(lik.nbSites_)
157  {
158  // Now adjust pointers:
159  transitionMatrix_->setHmmStateAlphabet(hiddenAlphabet_.get());
160  emissionProbabilities_->setHmmStateAlphabet(hiddenAlphabet_.get());
161  }
162 
164  {
166  AbstractParametrizable::operator=(lik);
167  hiddenAlphabet_ = std::unique_ptr<HmmStateAlphabet>(dynamic_cast<HmmStateAlphabet*>(lik.hiddenAlphabet_->clone()));
168  transitionMatrix_ = std::unique_ptr<HmmTransitionMatrix>(dynamic_cast<HmmTransitionMatrix*>(lik.transitionMatrix_->clone()));
169  emissionProbabilities_ = std::unique_ptr<HmmEmissionProbabilities>(dynamic_cast<HmmEmissionProbabilities*>(lik.emissionProbabilities_->clone()));
170  likelihood_ = lik.likelihood_;
175  scales_ = lik.scales_;
176  dScales_ = lik.dScales_;
177  d2Scales_ = lik.d2Scales_;
178  logLik_ = lik.logLik_;
180  nbStates_ = lik.nbStates_;
181  nbSites_ = lik.nbSites_;
182 
183  // Now adjust pointers:
184  transitionMatrix_->setHmmStateAlphabet(hiddenAlphabet_.get());
185  emissionProbabilities_->setHmmStateAlphabet(hiddenAlphabet_.get());
186  return *this;
187  }
188 
190 
191  RescaledHmmLikelihood* clone() const { return new RescaledHmmLikelihood(*this); }
192 
193 public:
196 
199 
202 
203  void setBreakPoints(const std::vector<size_t>& breakPoints)
204  {
205  breakPoints_ = breakPoints;
206  computeForward_();
207  backLikelihoodUpToDate_ = false;
208  }
209 
210  const std::vector<size_t>& getBreakPoints() const { return breakPoints_; }
211 
212  void setParameters(const ParameterList& pl)
213  {
215  }
216 
217  double getValue() const { return -logLik_; }
218 
219  double getLogLikelihood() const { return logLik_; }
220 
221  double getLikelihoodForASite(size_t site) const;
222 
223  double getDLogLikelihoodForASite(size_t site) const;
224 
225  double getD2LogLikelihoodForASite(size_t site) const;
226 
228 
229  void setNamespace(const std::string& nameSpace);
230 
231  void fireParameterChanged(const ParameterList& pl);
232 
234 
235  void getHiddenStatesPosteriorProbabilities(std::vector< std::vector<double> >& probs, bool append = false) const;
236 
237 protected:
238  void computeForward_();
239  void computeBackward_() const;
240 
241 
242  void computeDLikelihood_() const
243  {
245  }
246 
247  void computeD2Likelihood_() const
248  {
250  }
251 
252  void computeDForward_() const;
253 
254  void computeD2Forward_() const;
255 };
256 }
257 #endif // BPP_NUMERIC_HMM_RESCALEDHMMLIKELIHOOD_H
AbstractHmmLikelihood & operator=(const AbstractHmmLikelihood &adhlik)
A partial implementation of the Parametrizable interface.
void setParametersValues(const ParameterList &parameters)
Update the parameters from parameters.
Interface for computing emission probabilities in a Hidden Markov Model.
Hidden states alphabet.
Describe the transition probabilities between hidden states of a Hidden Markov Model.
The parameter list object.
Definition: ParameterList.h:65
A simple implementation of hidden Markov models recursion.
const HmmEmissionProbabilities & getHmmEmissionProbabilities() const
std::vector< double > dScales_
const HmmStateAlphabet & getHmmStateAlphabet() const
std::vector< double > likelihood_
The likelihood arrays.
void setParameters(const ParameterList &pl)
Set the point where the function must be computed.
std::vector< std::vector< double > > d2Likelihood_
std::vector< size_t > breakPoints_
const HmmTransitionMatrix & getHmmTransitionMatrix() const
std::vector< double > d2Scales_
double getDLogLikelihoodForASite(size_t site) const
void setBreakPoints(const std::vector< size_t > &breakPoints)
HmmStateAlphabet & getHmmStateAlphabet()
RescaledHmmLikelihood * clone() const
Create a copy of this object and send a pointer to it.
std::vector< std::vector< double > > dLikelihood_
derivatec of forward likelihood
Vdouble getHiddenStatesPosteriorProbabilitiesForASite(size_t site) const
double getValue() const
Get the value of the function at the current point.
RescaledHmmLikelihood(HmmStateAlphabet *hiddenAlphabet, HmmTransitionMatrix *transitionMatrix, HmmEmissionProbabilities *emissionProbabilities, const std::string &prefix)
Build a new RescaledHmmLikelihood object.
HmmEmissionProbabilities & getHmmEmissionProbabilities()
std::unique_ptr< HmmTransitionMatrix > transitionMatrix_
void getHiddenStatesPosteriorProbabilities(std::vector< std::vector< double > > &probs, bool append=false) const
std::unique_ptr< HmmEmissionProbabilities > emissionProbabilities_
double getD2LogLikelihoodForASite(size_t site) const
double getLikelihoodForASite(size_t site) const
Get the likelihood for a site, and its derivatives.
HmmTransitionMatrix & getHmmTransitionMatrix()
const std::vector< size_t > & getBreakPoints() const
RescaledHmmLikelihood & operator=(const RescaledHmmLikelihood &lik)
void fireParameterChanged(const ParameterList &pl)
Notify the class when one or several parameters have changed.
std::vector< double > scales_
scales for likelihood computing
void setNamespace(const std::string &nameSpace)
Set the namespace for the parameter names.
std::vector< std::vector< double > > backLikelihood_
backward likelihood
Vdouble getLikelihoodForEachSite() const
Get the likelihood for each site.
std::unique_ptr< HmmStateAlphabet > hiddenAlphabet_
The alphabet describing the hidden states.
RescaledHmmLikelihood(const RescaledHmmLikelihood &lik)
std::vector< double > Vdouble
Definition: VectorTools.h:70