bpp-phyl3  3.0.0
AlignedPhyloLikelihoodHmm.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_HMMOFALIGNEDPHYLOLIKELIHOOD_H
6 #define BPP_PHYL_LIKELIHOOD_PHYLOLIKELIHOODS_HMMOFALIGNEDPHYLOLIKELIHOOD_H
7 
8 
9 #include "HmmLikelihood_DF.h"
11 
12 // From Numeric
14 
15 #include "HmmPhyloAlphabet.h"
17 
18 
19 namespace bpp
20 {
32 {
33 private:
34  std::shared_ptr<HmmPhyloAlphabet> hma_;
35 
36  std::shared_ptr<FullHmmTransitionMatrix> htm_;
37 
38  std::shared_ptr<HmmPhyloEmissionProbabilities> hpep_;
39 
40  mutable std::shared_ptr<HmmLikelihood_DF> hmm_;
41 
42 public:
45  std::shared_ptr<PhyloLikelihoodContainer> pC,
46  const std::vector<size_t>& nPhylo,
47  bool inCollection = true);
48 
50 
51 protected:
58  hma_(mlc.hma_),
59  htm_(mlc.htm_),
60  hpep_(mlc.hpep_),
61  hmm_(mlc.hmm_)
62  {}
63 
65  {
67  hma_ = mlc.hma_;
68  htm_ = mlc.htm_;
69  hpep_ = mlc.hpep_;
70  hmm_ = mlc.hmm_;
71  return *this;
72  }
73 
74  AlignedPhyloLikelihoodHmm* clone() const override
75  {
76  return new AlignedPhyloLikelihoodHmm(*this);
77  }
78 
79 public:
80  void setNamespace(const std::string& nameSpace) override;
81 
82  void fireParameterChanged(const ParameterList& parameters) override;
83 
85  {
86  return *hma_;
87  }
88 
89  std::shared_ptr<const HmmPhyloAlphabet> getHmmStateAlphabet() const
90  {
91  return hma_;
92  }
93 
100  {
101  return *hmm_;
102  }
103 
104  std::shared_ptr<LikelihoodCalculation> getLikelihoodCalculation () const override
105  {
106  return hmm_;
107  }
108 
110  {
111  return *hmm_;
112  }
113 
114  std::shared_ptr<AlignedLikelihoodCalculation> getAlignedLikelihoodCalculation () const override
115  {
116  return hmm_;
117  }
118 
119  /*
120  *@brief return the posterior probabilities of subprocess on each site.
121  *
122  *@return MatrixXd sites x states
123  */
125  {
126  VVdouble pp;
127  auto mat = hmm_->getHiddenStatesPosteriorProbabilities().transpose();
128  copyEigenToBpp(mat, pp);
129  return pp;
130  }
131 
133  {
134  Vdouble pp;
135  auto mat = hmm_->getHiddenStatesPosteriorProbabilitiesForASite(site);
136  copyEigenToBpp(mat, pp);
137  return pp;
138  }
139 
140  const Eigen::MatrixXd& getHmmTransitionMatrix() const
141  {
142  return hmm_->getHmmTransitionMatrix();
143  }
144 };
145 } // end of namespace bpp.
146 #endif // BPP_PHYL_LIKELIHOOD_PHYLOLIKELIHOODS_HMMOFALIGNEDPHYLOLIKELIHOOD_H
The AlignedPhyloLikelihoodSet abstract class.
AbstractAlignedPhyloLikelihoodSet & operator=(const AbstractAlignedPhyloLikelihoodSet &soap)
The PhyloLikelihoodSet class, to manage a subset of PhyloLikelihoods from a given PhyloLikelihoodCont...
const Context & context() const override
Likelihood framework based on a hmm of simple likelihoods.
std::shared_ptr< HmmPhyloEmissionProbabilities > hpep_
AlignedLikelihoodCalculation & alignedLikelihoodCalculation() const override
std::shared_ptr< LikelihoodCalculation > getLikelihoodCalculation() const override
AlignedPhyloLikelihoodHmm(Context &context, std::shared_ptr< PhyloLikelihoodContainer > pC, const std::vector< size_t > &nPhylo, bool inCollection=true)
void fireParameterChanged(const ParameterList &parameters) override
Vdouble getPosteriorProbabilitiesForASitePerAligned(size_t site) const
AlignedPhyloLikelihoodHmm & operator=(const AlignedPhyloLikelihoodHmm &mlc)
std::shared_ptr< FullHmmTransitionMatrix > htm_
std::shared_ptr< const HmmPhyloAlphabet > getHmmStateAlphabet() const
const HmmPhyloAlphabet hmmStateAlphabet() const
std::shared_ptr< HmmLikelihood_DF > hmm_
VVdouble getPosteriorProbabilitiesPerSitePerAligned() const
std::shared_ptr< HmmPhyloAlphabet > hma_
void setNamespace(const std::string &nameSpace) override
const Eigen::MatrixXd & getHmmTransitionMatrix() const
std::shared_ptr< AlignedLikelihoodCalculation > getAlignedLikelihoodCalculation() const override
AlignedPhyloLikelihoodHmm * clone() const override
LikelihoodCalculation & likelihoodCalculation() const override
AlignedPhyloLikelihoodHmm(const AlignedPhyloLikelihoodHmm &mlc)
Context for dataflow node construction.
Definition: DataFlow.h:527
Hidden states alphabet.
Defines the basic types of data flow nodes.
std::vector< double > Vdouble
template void copyEigenToBpp(const ExtendedFloatMatrixXd &eigenMatrix, std::vector< std::vector< double >> &bppMatrix)
std::vector< Vdouble > VVdouble