bpp-phyl3 3.0.0
ForwardLikelihoodTree.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#include "Model.h"
7#include "Parametrizable.h"
8#include "Sequence_DF.h"
9
10
11#include <string> // debug
12
13
14using namespace bpp;
15using namespace std;
16
17ConditionalLikelihoodForwardRef ForwardLikelihoodTree::makeInitialConditionalLikelihood(
18 const string& sequenceName,
19 const AlignmentDataInterface& sites)
20{
21 size_t nbSites = sites.getNumberOfSites();
22
23 const auto sequenceIndex = sites.getSequencePosition (sequenceName);
24 Eigen::MatrixXd initCondLik ((int)nbState_, (int)nbSites);
25
26 for (size_t site = 0; site < nbSites; ++site)
27 {
28 for (auto state = 0; state < nbState_; ++state)
29 {
30 initCondLik (Eigen::Index (state), Eigen::Index (site)) =
31 sites (site, sequenceIndex, statemap_.getAlphabetStateAsInt(size_t(state)));
32 }
33 }
34
35 return Sequence_DF::create (context_, std::move(initCondLik), sequenceName);
36}
37
39 shared_ptr<ProcessEdge> processEdge,
40 const AlignmentDataInterface& sites)
41{
42 const auto brlen = processEdge->getBrLen();
43 const auto model = processEdge->getModel();
44 const auto nMod = processEdge->getNMod();
45 const auto brprob = processEdge->getProba();
46
47 auto childConditionalLikelihood = makeForwardLikelihoodAtNode (processTree_->getSon(processEdge), sites);
48
49 ForwardLikelihoodBelowRef forwardEdge;
50
51 auto zero = context_.getZero();
52
53 if (brlen) // Branch with transition through a model
54 {
55 if (dynamic_pointer_cast<const TransitionModelInterface>(model->targetValue()))
56 {
57 auto transitionMatrix = ConfiguredParametrizable::createMatrix<ConfiguredModel, TransitionMatrixFromModel, Eigen::MatrixXd>(context_, {model, brlen, zero, nMod}, transitionMatrixDimension (size_t(nbState_)));
58
59 processEdge->setTransitionMatrix(transitionMatrix);
60
61 forwardEdge = ForwardTransition::create (
62 context_, {transitionMatrix, childConditionalLikelihood}, likelihoodMatrixDim_);
63 }
64 else
65 {
66 auto transitionFunction = TransitionFunctionFromModel::create(context_, {model, brlen, zero}, transitionFunctionDimension(nbState_));
67
68 forwardEdge = ForwardTransitionFunction::create(context_, {childConditionalLikelihood, transitionFunction}, likelihoodMatrixDim_);
69 }
70 }
71 else if (brprob)
72 {
73 forwardEdge = ForwardProportion::create(
74 context_, {brprob, childConditionalLikelihood}, likelihoodMatrixDim_);
75 }
76 else // junction branch above a mixture node
77 {
78 forwardEdge = childConditionalLikelihood;
79 }
80
81 if (!hasEdgeIndex(forwardEdge)) // ie this edge does not exist before
82 {
83 // false (root) top-node of the edge (only way to build an edge in
84 // the DAG). Correct top-node will be set afterwards.
85
86 link(getRoot(), childConditionalLikelihood, forwardEdge);
87
88 setEdgeIndex(forwardEdge, processTree_->getEdgeIndex(processEdge)); // gets the index of the corresponding branch in the processTree_
89
90 auto spIndex = processEdge->getSpeciesIndex();
91
92 if (brprob || brlen)
93 {
94 if (mapEdgesIndexes_.find(spIndex) == mapEdgesIndexes_.end())
95 mapEdgesIndexes_[spIndex] = DAGindexes();
96 mapEdgesIndexes_[spIndex].push_back(getEdgeIndex(forwardEdge));
97 }
98 }
99
100#ifdef DEBUG
101 std::cerr << "E " << getEdgeIndex(forwardEdge) << " : forwardEdge " << forwardEdge << endl;
102 std::cerr << " -> N " << processTree_->getNodeIndex(processTree_->getSon(processEdge)) << std::endl;
103#endif
104
105 return forwardEdge;
106}
107
109{
110 const auto childBranches = processTree_->getBranches (processNode);
111
112 auto spIndex = processNode->getSpeciesIndex();
113
115
116 if (childBranches.empty ())
117 {
118 forwardNode = makeInitialConditionalLikelihood (processNode->getName (), sites);
119
120 if (!hasNodeIndex(forwardNode))
121 {
122 createNode(forwardNode);
123 setNodeIndex(forwardNode, processTree_->getNodeIndex(processNode));
124 if (mapNodesIndexes_.find(spIndex) == mapNodesIndexes_.end())
125 mapNodesIndexes_[spIndex] = DAGindexes();
126 mapNodesIndexes_[spIndex].push_back(getNodeIndex(forwardNode));
127 }
128 }
129 else
130 {
131 // depE are edges used to link ForwardLikelihoodTree edges to this
132 // node
133 std::vector<ForwardLikelihoodBelowRef> depE(childBranches.size());
134 NodeRefVec deps(childBranches.size());
135
136 for (size_t i = 0; i < childBranches.size (); ++i)
137 {
138 depE[i] = makeForwardLikelihoodAtEdge (childBranches[i], sites);
139 deps[i] = depE[i];
140 }
141
142 if (processNode->isSpeciation())
143 forwardNode = SpeciationForward::create(context_, std::move(deps),
145 else if (processNode->isMixture())
146 forwardNode = MixtureForward::create(context_, std::move(deps),
148 else
149 throw Exception("ForwardLikelihoodTree::makeConditionalLikelihoodAtNode : event not recognized for node " + TextTools::toString(processNode->getSpeciesIndex()));
150
151 // Fix the DAG
152 if (!hasNodeIndex(forwardNode))
153 {
154 createNode(forwardNode);
155 setNodeIndex(forwardNode, processTree_->getNodeIndex(processNode));
156
157 // Put the node in speciesIndex map if it is not son of a
158 // mixture node (which would hold the species index)
159
160 bool fathmixture(false);
161
162 if (processTree_->hasFather(processNode))
163 {
164 auto fatherNode = processTree_->getFatherOfNode (processNode);
165 fathmixture = fatherNode->isMixture();
166 }
167
168 if (!fathmixture)
169 {
170 if (mapNodesIndexes_.find(spIndex) == mapNodesIndexes_.end())
171 mapNodesIndexes_[spIndex] = DAGindexes();
172 mapNodesIndexes_[spIndex].push_back(getNodeIndex(forwardNode));
173 }
174
175 for (size_t i = 0; i < depE.size (); ++i)
176 {
177 auto fs = getNodes(depE[i]);
178 // fix the top node of the edge (because it was set to bidon
179 unlink(fs.first, fs.second);
180 link(forwardNode, fs.second, depE[i]);
181 }
182 }
183 }
184
185#ifdef DEBUG
186 std::cerr << "N " << getNodeIndex(forwardNode) << " : forwardNode " << forwardNode << endl;
187 for (size_t i = 0; i < childBranches.size (); ++i)
188 {
189 std::cerr << " -> E " << processTree_->getEdgeIndex(childBranches[i]) << std::endl;
190 }
191
192#endif
193
194 return forwardNode;
195}
196
197/************************************************************/
198/************************************************************
199 * For ProbabilityDAG
200 */
201
202ProbabilityDAG::ProbabilityDAG(std::shared_ptr<ForwardLikelihoodTree> forwardTree) :
203 DAProb(forwardTree->getGraph()), context_(forwardTree->getContext())
204{
206
207 associateNode(rootProb, forwardTree->getNodeGraphid(forwardTree->getRoot()));
208 setNodeIndex(rootProb, forwardTree->getRootIndex());
209
210 auto allIndex = forwardTree->getAllNodesIndexes();
211 std::vector<const Node_DF*> vN;
212
213 for (auto id: allIndex)
214 {
215 // Start from external nodes (which may be not leaves)
216 if (forwardTree->getOutgoingEdges(id).size() == 0)
217 {
218 auto n = makeProbaAtNode_(id, forwardTree).get();
219 n->targetValue();
220 vN.push_back(n);
221 }
222 }
223
224 // using bpp::DotOptions;
225 // writeGraphToDot("proba.dot", vN, DotOptions::DetailedNodeInfo);
226}
227
228
229ProbaRef ProbabilityDAG::makeProbaAtEdge_ (PhyloTree::EdgeIndex edgeIndex, std::shared_ptr<ForwardLikelihoodTree> forwardTree)
230{
231 auto fatherIndex = forwardTree->getFatherOfEdge(edgeIndex);
232 auto edgeForward = forwardTree->getEdge(edgeIndex);
233
234 // get/check if node with backward likelihood exists
235 auto probaNode = hasNode(fatherIndex)
236 ? getNode(fatherIndex)
237 : makeProbaAtNode_(fatherIndex, forwardTree);
238
239 const auto processEdge = forwardTree->getProcessTree()->getEdge(edgeIndex);
240
241 const auto brprob = processEdge->getProba();
242
243 ProbaRef probaEdge;
244
245 if (brprob)
246 {
247 probaEdge = ProbaMul::create(context_, {brprob, probaNode}, Dimension<double>());
248 }
249 else
250 {
251 probaEdge = Identity<double>::create(context_, {probaNode}, Dimension<double>(), edgeIndex);
252 }
253
254 // put object in the tree
255 if (!hasEdge(probaEdge))
256 {
257 associateEdge(probaEdge, forwardTree->getEdgeGraphid(edgeForward));
258 setEdgeIndex(probaEdge, edgeIndex);
259 }
260
261 return probaEdge;
262}
263
264
265ProbaRef ProbabilityDAG::makeProbaAtNode_ (PhyloTree::NodeIndex nodeIndex, std::shared_ptr<ForwardLikelihoodTree> forwardTree)
266{
267 // else get incoming edges
268 const auto edgesIndexes = forwardTree->getIncomingEdges(nodeIndex);
269
270 // get upper dependencies
271 NodeRefVec deps;
272
273 ProbaRef probaNode;
274 auto forwardNode = forwardTree->getNode(nodeIndex);
275
276 for (const auto& edgeIndex:edgesIndexes)
277 {
278 auto backEdge = hasEdge(edgeIndex)
279 ? getEdge(edgeIndex)
280 : makeProbaAtEdge_(edgeIndex, forwardTree);
281
282 if (edgesIndexes.size() == 1)
283 {
284 probaNode = Identity<double>::create(context_, {backEdge}, Dimension<double>(), nodeIndex);
285 break;
286 }
287 else
288 deps.push_back(backEdge);
289 }
290
291 // Then compute sum if several incoming nodes
292 if (deps.size() > 1)
293 probaNode = ProbaSum::create(context_, std::move(deps), Dimension<double>());
294
295 if (!hasNode(probaNode))
296 {
297 associateNode(probaNode, forwardTree->getNodeGraphid(forwardNode));
298 setNodeIndex(probaNode, nodeIndex);
299 }
300
301 return probaNode;
302}
AssociationGraphObserver< N, E >::EdgeIndex EdgeIndex
std::pair< Nref, Nref > getNodes(Eref edge) const
AssociationGraphObserver< N, E >::NodeIndex NodeIndex
virtual std::shared_ptr< N > getNode(NodeIndex nodeIndex) const=0
virtual void link(std::shared_ptr< N > nodeObjectA, std::shared_ptr< N > nodeObjectB, std::shared_ptr< E > edgeObject=00)=0
virtual NodeIndex setNodeIndex(const std::shared_ptr< N > nodeObject, NodeIndex index)=0
virtual void associateNode(std::shared_ptr< N > nodeObject, NodeGraphid node)=0
virtual std::shared_ptr< N > getRoot() const=0
virtual EdgeIndex getEdgeIndex(const std::shared_ptr< E > edgeObject) const=0
virtual bool hasNodeIndex(const std::shared_ptr< N > nodeObject) const=0
virtual void unlink(std::shared_ptr< N > nodeObjectA, std::shared_ptr< N > nodeObjectB)=0
virtual EdgeIndex setEdgeIndex(const std::shared_ptr< E > edgeObject, EdgeIndex index)=0
virtual NodeIndex getNodeIndex(const std::shared_ptr< N > nodeObject) const=0
virtual bool hasEdgeIndex(const std::shared_ptr< E > edgeObject) const=0
virtual bool hasNode(NodeIndex nodeIndex) const=0
virtual std::shared_ptr< E > getEdge(EdgeIndex edgeIndex) const=0
virtual bool hasEdge(EdgeIndex edgeIndex) const=0
virtual void associateEdge(std::shared_ptr< E > edgeObject, EdgeGraphid edge)=0
virtual void createNode(std::shared_ptr< N > newNodeObject)=0
static ValueRef< R > create(Context &c, NodeRefVec &&deps, const Dimension< R > &dim)
Build a new CWiseAdd node with the given output dimensions.
static ValueRef< R > create(Context &c, NodeRefVec &&deps, const Dimension< R > &dim)
Build a new CWiseApply node.
static std::shared_ptr< Self > create(Context &c, const Dimension< T > &dim)
Build a new ConstantOne node of the given dimension.
NodeRef getZero()
Definition: DataFlow.h:568
std::shared_ptr< ProcessTree > processTree_
ConditionalLikelihoodForwardRef makeInitialConditionalLikelihood(const std::string &sequenceName, const AlignmentDataInterface &sites)
Compute ConditionalLikelihood for leaf.
const StateMapInterface & statemap_
std::map< Speciesindex, DAGindexes > mapNodesIndexes_
ForwardLikelihoodBelowRef makeForwardLikelihoodAtEdge(std::shared_ptr< ProcessEdge > edge, const AlignmentDataInterface &sites)
Compute ConditionalLikelihood after reading edge on the forward process (ie at top of the edge).
std::map< Speciesindex, DAGindexes > mapEdgesIndexes_
ConditionalLikelihoodForwardRef makeForwardLikelihoodAtNode(std::shared_ptr< ProcessNode > node, const AlignmentDataInterface &sites)
static ValueRef< R > create(Context &c, NodeRefVec &&deps, const Dimension< R > &dim, size_t id=0)
Build a new Convert node with the given output dimensions.
static ValueRef< R > create(Context &c, NodeRefVec &&deps, const Dimension< R > &dim)
Build a new MatrixProduct node with the given output dimensions.
ProbaRef makeProbaAtNode_(PhyloTree::EdgeIndex edgeIndex, std::shared_ptr< ForwardLikelihoodTree > tree)
ProbabilityDAG(std::shared_ptr< ForwardLikelihoodTree > tree)
ProbaRef makeProbaAtEdge_(PhyloTree::EdgeIndex edgeIndex, std::shared_ptr< ForwardLikelihoodTree > tree)
computation of the probabilities with the same approach as for BackwardLikelihoodTree.
static ValueRef< T > create(Context &c, T &&value, const std::string &name)
Definition: Sequence_DF.h:31
virtual int getAlphabetStateAsInt(size_t index) const =0
virtual size_t getNumberOfSites() const=0
virtual size_t getSequencePosition(const std::string &sequenceKey) const =0
static std::shared_ptr< Self > create(Context &c, NodeRefVec &&deps, const Dimension< T > &dim)
Definition: Model.cpp:273
std::string toString(T t)
T zero(const Dimension< T > &)
Defines the basic types of data flow nodes.
Dimension< TransitionFunction > transitionFunctionDimension(Eigen::Index nbState)
std::vector< uint > DAGindexes
Helper: create a map with mutable dataflow nodes for each branch of the tree. The map is indexed by b...
ValueRef< MatrixLik > ConditionalLikelihoodForwardRef
std::vector< NodeRef > NodeRefVec
Alias for a dependency vector (of NodeRef).
Definition: DataFlow.h:81
std::shared_ptr< Value< double > > ProbaRef
ValueRef< MatrixLik > ForwardLikelihoodBelowRef
MatrixDimension transitionMatrixDimension(std::size_t nbState)
Definition: Model.h:19
Specialisation of Dimension<T> for floating point types.