bpp-phyl3 3.0.0
HmmProcessPhyloLikelihood.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_HMMPROCESSPHYLOLIKELIHOOD_H
6#define BPP_PHYL_LIKELIHOOD_PHYLOLIKELIHOODS_HMMPROCESSPHYLOLIKELIHOOD_H
7
8
9#include "../HmmSequenceEvolution.h"
12
13// From SeqLib:
15
16
17#include "HmmLikelihood_DF.h"
19
20
21namespace bpp
22{
32{
33private:
34 std::shared_ptr<HmmPhyloEmissionProbabilities> Hpep_;
35
39 mutable std::shared_ptr<HmmLikelihood_DF> hmm_;
40
41public:
43 std::shared_ptr<const AlignmentDataInterface> data,
44 std::shared_ptr<HmmSequenceEvolution> processSeqEvol,
45 std::shared_ptr<CollectionNodes> collNodes,
46 size_t nSeqEvol = 0,
47 size_t nData = 0);
48
49protected:
57 Hpep_(std::shared_ptr<HmmPhyloEmissionProbabilities>(mlc.Hpep_->clone())),
58 hmm_(std::shared_ptr<HmmLikelihood_DF>(mlc.hmm_->clone())) {}
59
60 HmmProcessPhyloLikelihood* clone() const override { return new HmmProcessPhyloLikelihood(*this); }
61
62public:
64
65public:
66 void setNamespace(const std::string& nameSpace) override;
67
68 void fireParameterChanged(const ParameterList& parameters) override;
69
76 {
77 return *hmm_;
78 }
79
80 std::shared_ptr<LikelihoodCalculation> getLikelihoodCalculation () const override
81 {
82 return hmm_;
83 }
84
86 {
87 return *hmm_;
88 }
89
90 std::shared_ptr<AlignedLikelihoodCalculation> getAlignedLikelihoodCalculation () const override
91 {
92 return hmm_;
93 }
94
95 /*
96 *@brief return the posterior probabilities of subprocess on each site.
97 *
98 *@return MatrixXd sites x states
99 */
101 {
102 VVdouble pp;
103 auto mat = hmm_->getHiddenStatesPosteriorProbabilities().transpose();
104 copyEigenToBpp(mat, pp);
105 return pp;
106 }
107
108 const Eigen::MatrixXd& getHmmTransitionMatrix() const
109 {
110 return hmm_->getHmmTransitionMatrix();
111 }
112};
113} // end of namespace bpp.
114#endif // BPP_PHYL_LIKELIHOOD_PHYLOLIKELIHOODS_HMMPROCESSPHYLOLIKELIHOOD_H
A simple implementation of hidden Markov models recursion, in DataFlow implementation.
Emission probabilities in the context of DF phylolikeihoods.
Likelihood framework based on a hmm of simple likelihoods.
std::shared_ptr< LikelihoodCalculation > getLikelihoodCalculation() const override
const Eigen::MatrixXd & getHmmTransitionMatrix() const
AlignedLikelihoodCalculation & alignedLikelihoodCalculation() const override
void fireParameterChanged(const ParameterList &parameters) override
HmmProcessPhyloLikelihood(std::shared_ptr< const AlignmentDataInterface > data, std::shared_ptr< HmmSequenceEvolution > processSeqEvol, std::shared_ptr< CollectionNodes > collNodes, size_t nSeqEvol=0, size_t nData=0)
std::shared_ptr< AlignedLikelihoodCalculation > getAlignedLikelihoodCalculation() const override
std::shared_ptr< HmmPhyloEmissionProbabilities > Hpep_
LikelihoodCalculation & likelihoodCalculation() const override
HmmProcessPhyloLikelihood(const HmmProcessPhyloLikelihood &mlc)
std::shared_ptr< HmmLikelihood_DF > hmm_
LikelihoodCalculation in context of HMM.
VVdouble getPosteriorProbabilitiesPerSitePerProcess() const override
HmmProcessPhyloLikelihood * clone() const override
void setNamespace(const std::string &nameSpace) override
Partial implementation of the Likelihood interface for multiple processes.
Defines the basic types of data flow nodes.
template void copyEigenToBpp(const ExtendedFloatMatrixXd &eigenMatrix, std::vector< std::vector< double > > &bppMatrix)
std::vector< Vdouble > VVdouble