bpp-phyl3 3.0.0
HmmLikelihoodComputation.cpp
Go to the documentation of this file.
1// SPDX-FileCopyrightText: The Bio++ Development Group
2//
3// SPDX-License-Identifier: CECILL-2.1
4
5#include "HmmLikelihood_DF.h"
6
7// from the STL:
8#include <iostream>
9#include <algorithm>
10using namespace bpp;
11using namespace std;
12using namespace numeric;
13
14/***********************************/
15/* Forward Likelihood */
16/***********************************/
17
18void ForwardHmmLikelihood_DF::compute()
19{
20 const auto& hmmEq = accessValueConstCast<Eigen::VectorXd>(*this->dependency(0));
21
22 const auto& hmmTrans = accessValueConstCast<Eigen::MatrixXd>(*this->dependency(1));
23 const auto& hmmEmis = accessValueConstCast<MatrixLik>(*this->dependency(2));
24
25 auto& condLik = dynamic_pointer_cast<CondLikelihood>(condLik_)->accessValueMutable();
26
27 auto nbSites = hmmEmis.cols();
28 auto nbStates = hmmEmis.rows();
29
30 VDataLik tscales((size_t)nbSites);
31
32 VectorLik tmp((int)nbStates);
33
34 // Initialisation:
35 parCondLik_[0] = hmmTrans.transpose() * hmmEq;
36
37 cwise(tmp) = cwise(parCondLik_[0]) * cwise(hmmEmis.col(0));
38 tscales[0] = tmp.sum();
39
40 for (auto s = 0; s < nbStates; s++)
41 {
42 condLik(s, 0) = convert(tmp(s) / tscales[0]);
43 }
44
45 // Iteration
46 for (auto i = 1; i < nbSites; i++)
47 {
48 parCondLik_[(size_t)i] = hmmTrans.transpose() * condLik.col(i - 1);
49
50 cwise(tmp) = cwise(hmmEmis.col(i)) * cwise(parCondLik_[(size_t)i]);
51 tscales[(size_t)i] = tmp.sum();
52
53 // tmp = condLik * scales
54 for (auto s = 0; s < nbStates; s++)
55 {
56 condLik(s, i) = convert(tmp(s) / tscales[(size_t)i]);
57 }
58 }
59 copyBppToEigen(tscales, this->accessValueMutable ());
60}
61
63{
64 if (&node == this)
65 {
67 }
68
69 /*
70 * 1st order derivatives of Forward Likelihood Arrays
71 *
72 * Dependencies are:
73 * Value<VectorXd> : Starting vector of states probabililies
74 * Value<MatrixXd> : TransitionMatrix
75 * Value<MatrixLik> : Matrix of Emission likelihoods states X sites
76 *
77 * ForwardHmmLikelihood_DF : Forward Computations
78 *
79 * Value<VectorXd> : Derivatives of starting vector of states probabililies
80 * Value<MatrixXd> : Derivatives of TransitionMatrix
81 * Value<MatrixLik> : Derivatives Matrix of Emission likelihoods states X sites
82 */
83
85 this->dependency(0),
86 this->dependency(1),
87 this->dependency(2),
88 this->shared_from_this(),
89 this->dependency(0)->derive (c, node),
90 this->dependency(1)->derive (c, node),
91 this->dependency(2)->derive (c, node)}, targetDimension_);
92}
93
94/***********************************/
95/* 1st order derivative of Forward Likelihood */
96/***********************************/
97
98
100{
101 const auto& hmmEq = accessValueConstCast<Eigen::VectorXd>(*this->dependency(0));
102
103 const auto& hmmTrans = accessValueConstCast<Eigen::MatrixXd>(*this->dependency(1));
104 const auto& hmmEmis = accessValueConstCast<MatrixLik>(*this->dependency(2));
105
106
107 auto forwardNode = dynamic_pointer_cast<ForwardHmmLikelihood_DF>(this->dependency(3));
108
109 const auto& condLik = forwardNode->getForwardCondLikelihood()->targetValue();
110
111 const auto& parCondLik = forwardNode->getParCondLik();
112
113 const auto& scales = forwardNode->targetValue();
114
115 const auto& dHmmEq = accessValueConstCast<Eigen::VectorXd>(*this->dependency(4));
116
117 const auto& dHmmTrans = accessValueConstCast<Eigen::MatrixXd>(*this->dependency(5));
118
119 const auto& dHmmEmis = accessValueConstCast<MatrixLik>(*this->dependency(6));
120
121 auto nbSites = hmmEmis.cols();
122 const int nbStates = (int)hmmEmis.rows();
123
124 VDataLik tdScales((size_t)nbSites);
125
126 auto& dCondLik = dynamic_pointer_cast<CondLikelihood>(dCondLik_)->accessValueMutable();
127
128 VectorLik dtmp(nbStates);
129
130 // Initialisation
131 dParCondLik_[0] = dHmmTrans.transpose() * hmmEq + hmmTrans.transpose() * dHmmEq;
132
133 cwise(dtmp) = (cwise(dParCondLik_[0]) * cwise(hmmEmis.col(0))
134 + cwise(dHmmEmis.col(0)) * cwise(parCondLik[0]));
135 tdScales[0] = dtmp.sum();
136
137 // dtmp = dCondLik * scales + CondLik * dScales
138
139 for (auto s = 0; s < nbStates; s++)
140 {
141 dCondLik(s, 0) = convert((dtmp(s) - condLik(s, 0) * tdScales[0]) / scales(0));
142 }
143
144 // Iteration
145 for (auto i = 1; i < nbSites; i++)
146 {
147 dParCondLik_[(size_t)i] = dHmmTrans.transpose() * condLik.col(i - 1) + hmmTrans.transpose() * dCondLik.col(i - 1);
148
149 cwise(dtmp) = cwise(dParCondLik_[(size_t)i]) * cwise(hmmEmis.col(i)) + cwise(parCondLik[(size_t)i]) * cwise(dHmmEmis.col(i));
150 tdScales[(size_t)i] = dtmp.sum();
151
152 for (auto s = 0; s < nbStates; s++)
153 {
154 dCondLik(s, i) = convert((dtmp(s) - condLik(s, i) * tdScales[(size_t)i]) / scales(i));
155 }
156 }
157
158 copyBppToEigen(tdScales, this->accessValueMutable ());
159}
160
162{
163 if (&node == this)
164 {
166 }
167
168 /*
169 *
170 * 2nd order derivatives of Forward Likelihood Arrays
171 *
172 * Dependencies are:
173 * Value<VectorXd> : Starting vector of states probabililies
174 * Value<MatrixXd> : TransitionMatrix
175 * Value<MatrixLik> : Matrix of Emission likelihoods states X sites
176 *
177 * ForwardHmmLikelihood_DF : Forward Computations
178 *
179 * Value<VectorXd> : 1st Derivatives of starting vector of states probabililies
180 * Value<MatrixXd> : 1st Derivatives of TransitionMatrix
181 * Value<MatrixLik> : 1st Derivatives Matrix of Emission likelihoods states X sites
182 *
183 * ForwardHmmDLikelihood_DF : 1st order derivatives Forward Computations
184 *
185 * Value<VectorXd> : 2nd Derivatives of starting vector of states probabililies
186 * Value<MatrixXd> : 2nd Derivatives of TransitionMatrix
187 * Value<MatrixLik> : 2nd Derivatives Matrix of Emission likelihoods states X sites
188 */
189
191 this->dependency(0),
192 this->dependency(1),
193 this->dependency(2),
194
195 this->dependency(3),
196
197 this->dependency(4),
198 this->dependency(5),
199 this->dependency(6),
200
201 this->shared_from_this(),
202
203 this->dependency(4)->derive (c, node),
204 this->dependency(5)->derive (c, node),
205 this->dependency(6)->derive (c, node)},
206
208}
209
210
211/***********************************/
212/* 2nd order derivative of Forward Likelihood */
213/***********************************/
214
215
217{
218 const auto& hmmEq = accessValueConstCast<Eigen::VectorXd>(*this->dependency(0));
219
220 const auto& hmmTrans = accessValueConstCast<Eigen::MatrixXd>(*this->dependency(1));
221 const auto& hmmEmis = accessValueConstCast<MatrixLik>(*this->dependency(2));
222
223
224 auto forwardNode = dynamic_pointer_cast<ForwardHmmLikelihood_DF>(this->dependency(3));
225
226 const auto& condLik = forwardNode->getForwardCondLikelihood()->targetValue();
227
228 const auto& parCondLik = forwardNode->getParCondLik();
229
230 const auto& scales = forwardNode->targetValue();
231
232
233 const auto& dHmmEq = accessValueConstCast<Eigen::VectorXd>(*this->dependency(4));
234
235 const auto& dHmmTrans = accessValueConstCast<Eigen::MatrixXd>(*this->dependency(5));
236
237 const auto& dHmmEmis = accessValueConstCast<MatrixLik>(*this->dependency(6));
238
239
240 auto forwardDNode = dynamic_pointer_cast<ForwardHmmDLikelihood_DF>(this->dependency(7));
241
242 const auto& dCondLik = forwardDNode->getForwardDCondLikelihood()->targetValue();
243
244 const auto& parDCondLik = forwardDNode->getParDCondLik();
245
246 const auto& dScales = forwardDNode->targetValue();
247
248
249 const auto& d2HmmEq = accessValueConstCast<Eigen::VectorXd>(*this->dependency(8));
250
251 const auto& d2HmmTrans = accessValueConstCast<Eigen::MatrixXd>(*this->dependency(9));
252
253 const auto& d2HmmEmis = accessValueConstCast<MatrixLik>(*this->dependency(10));
254
255
257
258 auto nbSites = hmmEmis.cols();
259 const int nbStates = static_cast<int>(hmmEmis.rows());
260
261 VDataLik td2Scales(static_cast<size_t>(nbSites));
262
263 Eigen::VectorXd d2CondLik(static_cast<int>(hmmEmis.rows()));
264
265 Eigen::VectorXd d2ParCondLik;
266 VectorLik d2tmp(nbStates);
267
268 // Initialisation:
269 d2ParCondLik = d2HmmTrans.transpose() * hmmEq + 2 * dHmmTrans.transpose() * dHmmEq + hmmTrans.transpose() * d2HmmEq;
270
271 cwise(d2tmp) = cwise(d2ParCondLik) * cwise(hmmEmis.col(0))
272 + 2 * cwise(parDCondLik[0]) * cwise(dHmmEmis.col(0))
273 + cwise(parCondLik[0]) * cwise(d2HmmEmis.col(0));
274
275 td2Scales[0] = d2tmp.sum();
276
277 // d2tmp = d2CondLik * scales + 2 * dCondLik * dScales + CondLik * d2Scales
278
279 for (auto s = 0; s < nbStates; s++)
280 {
281 d2CondLik(s, 0) = convert((d2tmp(s) - 2 * dCondLik(s, 0) * dScales(0) - condLik(s, 0) * td2Scales[0]) / scales(0));
282 }
283
284 // Iteration
285 for (auto i = 1; i < nbSites; i++)
286 {
287 d2ParCondLik = d2HmmTrans.transpose() * condLik.col(i)
288 + 2 * dHmmTrans.transpose() * dCondLik.col(i) + hmmTrans.transpose() * d2CondLik;
289
290 cwise(d2tmp) = cwise(d2ParCondLik) * cwise(hmmEmis.col(i))
291 + 2 * cwise(parDCondLik[(size_t)i]) * cwise(dHmmEmis.col(i))
292 + cwise(parCondLik[(size_t)i]) * cwise(d2HmmEmis.col(i));
293
294 td2Scales[(size_t)i] = d2tmp.sum();
295
296 for (auto s = 0; s < nbStates; s++)
297 {
298 d2CondLik(s, i) = convert((d2tmp(s) - 2 * dCondLik(s, i) * dScales(i) - condLik(s, i) * td2Scales[(size_t)i]) / scales(i));
299 }
300 }
301
302 copyBppToEigen(td2Scales, this->accessValueMutable ());
303}
304
305/*****************************************
306** Backward
307*****************************************/
309{
310 const auto& scales = accessValueConstCast<RowLik>(*this->dependency(0));
311
312 const auto& hmmTrans = accessValueConstCast<Eigen::MatrixXd>(*this->dependency(1));
313
314 const auto& hmmEmis = accessValueConstCast<MatrixLik>(*this->dependency(2));
315
316 auto nbSites = hmmEmis.cols();
317 auto nbStates = hmmEmis.rows();
318
319 std::vector<Eigen::VectorXd> tScales((size_t)nbSites);
320
321 VectorLik tmp((int)nbStates);
322
323 // Initialisation:
324 tScales[(size_t)(nbSites - 1)] = Eigen::VectorXd::Constant(nbStates, 1);
325
326 // Iteration
327 for (auto i = nbSites - 1; i > 0; i--)
328 {
329 tScales[(size_t)(i - 1)].resize(nbStates);
330
331 cwise(tmp) = cwise(hmmEmis.col(i)) * cwise(tScales[(size_t)i]);
332
333 auto tmp2 = hmmTrans * tmp;
334
335 for (auto s = 0; s < nbStates; s++)
336 {
337 tScales[(size_t)(i - 1)](s) = convert(tmp2(s) / scales(i));
338 }
339 }
340
341 copyBppToEigen(tScales, this->accessValueMutable ());
342}
void compute() override
Computation implementation.
static std::shared_ptr< Self > create(Context &c, const Dimension< T > &dim)
Build a new ConstantOne node of the given dimension.
Context for dataflow node construction.
Definition: DataFlow.h:527
const ExtendedFloat & sum() const
static ValueRef< RowLik > create(Context &c, NodeRefVec &&deps, const Dimension< Eigen::MatrixXd > &dim)
void compute() override
Computation implementation.
Dimension< Eigen::MatrixXd > targetDimension_
Dimension of the data : states X sites.
std::vector< Eigen::VectorXd > dParCondLik_
static ValueRef< RowLik > create(Context &c, NodeRefVec &&deps, const Dimension< Eigen::MatrixXd > &dim)
void compute() override
Computation implementation.
NodeRef derive(Context &c, const Node_DF &node) final
Returns a node computing d(this_node_expression)/d(node_expression).
ValueRef< Eigen::MatrixXd > dCondLik_
derivatives of the conditional forward likelihoods : Will be used by likelihoods computation of 2nd o...
Dimension< Eigen::MatrixXd > targetDimension_
std::vector< Eigen::VectorXd > parCondLik_
ValueRef< Eigen::MatrixXd > condLik_
NodeRef derive(Context &c, const Node_DF &node) final
Returns a node computing d(this_node_expression)/d(node_expression).
Base dataflow Node class.
Definition: DataFlow.h:152
const NodeRef & dependency(std::size_t i) const noexcept
Definition: DataFlow.h:185
RowLik & accessValueMutable() noexcept
Definition: DataFlow.h:416
T & cwise(T &t)
Definition: DataFlowCWise.h:29
Defines the basic types of data flow nodes.
std::vector< DataLik > VDataLik
Definition: Definitions.h:23
template void copyBppToEigen(const std::vector< ExtendedFloat > &bppVector, Eigen::RowVectorXd &eigenVector)
double convert(const bpp::ExtendedFloat &ef)
std::shared_ptr< Node_DF > NodeRef
Definition: DataFlow.h:78