bpp-phyl3  3.0.0
AutoCorrelationProcessPhyloLikelihood.cpp
Go to the documentation of this file.
1 // SPDX-FileCopyrightText: The Bio++ Development Group
2 //
3 // SPDX-License-Identifier: CECILL-2.1
4 
6 
8 #include "HmmLikelihood_DF.h"
9 
10 using namespace std;
11 using namespace bpp;
12 
13 /******************************************************************************/
14 
15 AutoCorrelationProcessPhyloLikelihood::AutoCorrelationProcessPhyloLikelihood(
16  std::shared_ptr<const AlignmentDataInterface> data,
17  std::shared_ptr<AutoCorrelationSequenceEvolution> processSeqEvol,
18  std::shared_ptr<CollectionNodes> collNodes,
19  size_t nSeqEvol,
20  size_t nData) :
21  AbstractPhyloLikelihood(collNodes->context()),
22  AbstractAlignedPhyloLikelihood(collNodes->context(), data->getNumberOfSites()),
24  collNodes->context(),
25  data->getNumberOfSites(),
26  (processSeqEvol->getSubstitutionProcessNumbers().size() != 0)
27  ? processSeqEvol->substitutionProcess(processSeqEvol->getSubstitutionProcessNumbers()[0]).getNumberOfStates()
28  : 0,
29  nData),
30  AbstractSequencePhyloLikelihood(collNodes->context(), processSeqEvol, nData),
32  MultiProcessSequencePhyloLikelihood(data, processSeqEvol, collNodes, nSeqEvol, nData),
33  Hpep_(),
34  hmm_()
35 {
36  auto alphyl = make_shared<HmmPhyloAlphabet>(*this);
37 
38  Hpep_ = make_shared<HmmPhyloEmissionProbabilities>(alphyl);
39 
40  hmm_ = shared_ptr<HmmLikelihood_DF>(new HmmLikelihood_DF(context(), processSeqEvol->getHmmProcessAlphabet(), processSeqEvol->getHmmTransitionMatrix(), Hpep_));
41 
42  addParameters_(hmm_->getParameters());
43 }
44 
45 void AutoCorrelationProcessPhyloLikelihood::setNamespace(const std::string& nameSpace)
46 {
48  hmm_->setNamespace(nameSpace);
49 }
50 
51 
53 {
55 
56  hmm_->matchParametersValues(parameters);
57 }
virtual void fireParameterChanged(const ParameterList &parameters) override
void setNamespace(const std::string &nameSpace) override
virtual void addParameters_(const ParameterList &parameters)
const Context & context() const override
void fireParameterChanged(const ParameterList &parameters) override
void setNamespace(const std::string &nameSpace) override
std::shared_ptr< HmmPhyloEmissionProbabilities > Hpep_
A simple implementation of hidden Markov models recursion, in DataFlow implementation.
Partial implementation of the Likelihood interface for multiple processes.
Defines the basic types of data flow nodes.