bpp-phyl3 3.0.0
HmmPhyloEmissionProbabilities.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: The Bio++ Development Group
2//
3// SPDX-License-Identifier: CECILL-2.1
4
5#ifndef BPP_PHYL_LIKELIHOOD_PHYLOLIKELIHOODS_HMMPHYLOEMISSIONPROBABILITIES_H
6#define BPP_PHYL_LIKELIHOOD_PHYLOLIKELIHOODS_HMMPHYLOEMISSIONPROBABILITIES_H
7
9
10#include "../DataFlow/DataFlowCWise.h"
11#include "../HmmEmissionProbabilities_Eigen.h"
12#include "HmmPhyloAlphabet.h"
13
14namespace bpp
15{
17
24 public virtual HmmEmissionProbabilities_Eigen,
26{
27private:
29
30 std::shared_ptr<HmmPhyloAlphabet> phylAlph_;
31
32 /*
33 *@brief Emission likelihoods are stored in a Matrix from a set of
34 * RowVectors.
35 *
36 */
37
39
40 size_t nbSites_;
41
42public:
43 HmmPhyloEmissionProbabilities(std::shared_ptr<HmmPhyloAlphabet> alphabet);
44
47 context_(hEP.context_),
49 emProb_(hEP.emProb_),
51 {}
52
54
56 {
57 return phylAlph_.get();
58 }
59
60 size_t getNumberOfStates() const
61 {
62 return phylAlph_->getNumberOfStates();
63 }
64
65 size_t getNumberOfSites() const
66 {
67 return nbSites_;
68 }
69
76 void setHmmStateAlphabet(std::shared_ptr<HmmStateAlphabet> stateAlphabet);
77
88 DataLik operator()(size_t pos, size_t state) const
89 {
90 return (emProb_->targetValue())(Eigen::Index(state), Eigen::Index(pos));
91 }
92
94 {
95 return emProb_;
96 }
97
107 VectorLik operator()(size_t pos) const
108 {
109 return emProb_->targetValue().col(Eigen::Index(pos));
110 }
111
115 size_t getNumberOfPositions() const
116 {
117 return nbSites_;
118 }
119};
120} // end of namespace bpp.
121#endif // BPP_PHYL_LIKELIHOOD_PHYLOLIKELIHOODS_HMMPHYLOEMISSIONPROBABILITIES_H
Context for dataflow node construction.
Definition: DataFlow.h:527
Interface for computing emission probabilities in a Hidden Markov Model.
Emission probabilities in the context of DF phylolikeihoods.
VectorLik operator()(size_t pos) const
Operator access to the emission probabilities.
DataLik operator()(size_t pos, size_t state) const
Operator access to the emission probabilities.
HmmPhyloEmissionProbabilities * clone() const
HmmPhyloEmissionProbabilities(std::shared_ptr< HmmPhyloAlphabet > alphabet)
HmmPhyloEmissionProbabilities(const HmmPhyloEmissionProbabilities &hEP)
void setHmmStateAlphabet(std::shared_ptr< HmmStateAlphabet > stateAlphabet)
Set the new hidden state alphabet.
const HmmStateAlphabet * getHmmStateAlphabet() const
std::shared_ptr< HmmPhyloAlphabet > phylAlph_
Defines the basic types of data flow nodes.
std::shared_ptr< Value< T > > ValueRef
Shared pointer alias for Value<T>.
Definition: DataFlow.h:84