bpp-phyl3 3.0.0
DataFlowCWiseComputing.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_DATAFLOW_DATAFLOWCWISECOMPUTING_H
6#define BPP_PHYL_LIKELIHOOD_DATAFLOW_DATAFLOWCWISECOMPUTING_H
7
11#include <Eigen/Core>
12#include <algorithm>
13#include <cassert>
14#include <iostream>
15#include <list>
16#include <string>
17#include <tuple>
18#include <type_traits>
19
20#include "DataFlowCWise.h"
21#include "DataFlowNumeric.h"
22#include "Definitions.h"
23
24namespace bpp
25{
26/******************************************************************************
27 * Data flow nodes for those numerical functions.
28 *
29 * TODO numerical simplification:
30 * add(x,x) -> 2*x ? (and similar for mul, ...)
31 * all deps constant => return constant ?
32 */
33
34template<typename Result, typename From> class CWiseAdd;
35template<typename Result, typename From, typename Prop> class CWiseMean;
36template<typename Result, typename From> class CWiseSub;
37template<typename Result, typename From> class CWiseMul;
38template<typename Result, typename From> class CWiseDiv;
39
40template<typename R, typename T0, typename T1> class MatrixProduct;
41// Return same type as given
42template<typename T> class CWiseNegate;
43template<typename T> class CWiseInverse;
44template<typename T> class CWiseLog;
45template<typename T> class CWiseExp;
46template<typename T> class CWiseConstantPow;
47
48// Return DataLik (R)
49template<typename R, typename T0, typename T1> class ScalarProduct;
50template<typename R, typename T0, typename T1> class LogSumExp;
51
52// Return double
53template<typename F> class SumOfLogarithms;
54
55// Return same type as given
56
57template<typename T> class ShiftDelta;
58template<typename T> class CombineDeltaShifted;
59
60
61/*************************************************************************
62 * @brief r = f(x0) for each component or column
63 * - r: R.
64 * - x0: R
65 * - f : function between columns or components of R (type F)
66 *
67 */
68
69template<typename R, typename T, typename F> class CWiseApply : public Value<R>
70{
71public:
73
75 static ValueRef<R> create (Context& c, NodeRefVec&& deps, const Dimension<R>& dim)
76 {
77 // Check dependencies
78 checkDependenciesNotNull (typeid (Self), deps);
79 checkDependencyVectorSize (typeid (Self), deps, 2);
80 checkNthDependencyIsValue<T>(typeid (Self), deps, 0);
81 checkNthDependencyIsValue<F>(typeid (Self), deps, 1);
82
83 return cachedAs<Value<R>>(c, std::make_shared<Self>(std::move (deps), dim));
84 }
85
86 CWiseApply (NodeRefVec&& deps, const Dimension<R>& dim)
87 : Value<R>(std::move (deps)), bppLik_(size_t(dim.cols)), targetDimension_ (dim)
88 {
89 this->accessValueMutable().resize(dim.rows, dim.cols);
90 }
91
92 std::string debugInfo () const override
93 {
94 using namespace numeric;
95 return debug (this->accessValueConst ()) + " targetDim=" + to_string (targetDimension_);
96 }
97
98 // CWiseApply additional arguments = ().
99 bool compareAdditionalArguments (const Node_DF& other) const final
100 {
101 return dynamic_cast<const Self*>(&other) != nullptr;
102 }
103
104 // df(X(theta))/dtheta|X(theta) = df/dtheta|X(theta) + df/dX.dX/dtheta|X(theta)
105 NodeRef derive (Context& c, const Node_DF& node) final
106 {
107 if (&node == this)
108 {
110 }
111
112 NodeRefVec dep(2);
113 NodeRef df = this->dependency(1)->derive (c, node);
114
115 if (df->hasNumericalProperty (NumericalProperty::ConstantZero))
117 else
118 dep[0] = Self::create (c, {this->dependency(0), df}, targetDimension_);
119
120 NodeRef dX = this->dependency(0)->derive (c, node);
121
122 if (dX->hasNumericalProperty (NumericalProperty::ConstantZero))
124 else // product df/dX.dX/dtheta|X(theta)
125 {
126 NodeRef dfX = this->dependency(1)->derive(c, NodeX);
127 NodeRef dfXX = Self::create(c, {this->dependency(0), dfX}, targetDimension_);
129 }
130
131 return CWiseAdd<R, std::tuple<R, R>>::create (c, {std::move(dep)}, targetDimension_);
132 }
133
135 {
136 return Self::create (c, std::move (deps), targetDimension_);
137 }
138
139 std::string shape() const override
140 {
141 return "doubleoctagon";
142 }
143
144 std::string color() const override
145 {
146 return "#5e5e5e";
147 }
148
149 std::string description() const override
150 {
151 return "Function Apply";
152 }
153
154private:
155 void compute() override { compute<T>();}
156
157 template<class U>
158 typename std::enable_if<std::is_base_of<U, MatrixLik>::value && std::is_same<F, TransitionFunction>::value, void>::type
160 {
161 using namespace numeric;
162 auto& result = this->accessValueMutable ();
163 const auto& x0 = accessValueConstCast<T>(*this->dependency (0));
164 const auto& func = accessValueConstCast<F>(*this->dependency (1));
165
166 for (auto i = 0; i < x0.cols(); i++)
167 {
168 bppLik_[(size_t)i] = func(x0.col(i));
169 }
170
171 copyBppToEigen(bppLik_, result);
172
173#ifdef DEBUG
174 std::cerr << "=== Function Apply === " << this << std::endl;
175 std::cerr << "x0= " << x0 << std::endl;
176 std::cerr << "res=" << result << std::endl;
177 std::cerr << "=== end Function Apply === " << this << std::endl << std::endl;
178#endif
179 }
180
186 std::vector<VectorLik> bppLik_;
187
189};
190
191
192/*************************************************************************
193 * @brief r = x0 + x1 for each component.
194 * - r: R.
195 * - x0: T0.
196 * - x1: T1.
197 *
198 * Values converted to R with the semantics of numeric::convert.
199 * Node construction should be done with the create static method.
200 *
201 * Only defined for N = 2 for now.
202 * The generic version is horrible in C++11 (lack of auto return).
203 * Generic simplification routine is horrible too.
204 */
205
206template<typename R, typename T0, typename T1> class CWiseAdd<R, std::tuple<T0, T1>> : public Value<R>
207{
208public:
209 using Self = CWiseAdd;
210
212 static ValueRef<R> create (Context& c, NodeRefVec&& deps, const Dimension<R>& dim)
213 {
214 // Check dependencies
215 checkDependenciesNotNull (typeid (Self), deps);
216 checkDependencyVectorSize (typeid (Self), deps, 2);
217 checkNthDependencyIsValue<T0>(typeid (Self), deps, 0);
218 checkNthDependencyIsValue<T1>(typeid (Self), deps, 1);
219 // Select node implementation
220 bool zeroDep0 = deps[0]->hasNumericalProperty (NumericalProperty::ConstantZero);
221 bool zeroDep1 = deps[1]->hasNumericalProperty (NumericalProperty::ConstantZero);
222 if (zeroDep0 && zeroDep1)
223 {
224 return ConstantZero<R>::create (c, dim);
225 }
226 else if (zeroDep0 && !zeroDep1)
227 {
228 return Convert<R, T1>::create (c, {deps[1]}, dim);
229 }
230 else if (!zeroDep0 && zeroDep1)
231 {
232 return Convert<R, T0>::create (c, {deps[0]}, dim);
233 }
234 else
235 {
236 return cachedAs<Value<R>>(c, std::make_shared<Self>(std::move (deps), dim));
237 }
238 }
239
240 CWiseAdd (NodeRefVec&& deps, const Dimension<R>& dim)
241 : Value<R>(std::move (deps)), targetDimension_ (dim) {}
242
243 std::string debugInfo () const override
244 {
245 using namespace numeric;
246 return debug (this->accessValueConst ()) + " targetDim=" + to_string (targetDimension_);
247 }
248
249 std::string shape() const override
250 {
251 return "triangle";
252 }
253
254 std::string color() const override
255 {
256 return "#ff6060";
257 }
258
259 std::string description() const override
260 {
261 return "Add";
262 }
263
264 // Convert<T> additional arguments = ().
265 bool compareAdditionalArguments (const Node_DF& other) const final
266 {
267 return dynamic_cast<const Self*>(&other) != nullptr;
268 }
269
270 NodeRef derive (Context& c, const Node_DF& node) final
271 {
272 if (&node == this)
273 {
275 }
276 constexpr std::size_t n = 2;
277 NodeRefVec derivedDeps (n);
278 for (std::size_t i = 0; i < n; ++i)
279 {
280 derivedDeps[i] = this->dependency (i)->derive (c, node);
281 }
282 return Self::create (c, std::move (derivedDeps), targetDimension_);
283 }
284
286 {
287 return Self::create (c, std::move (deps), targetDimension_);
288 }
289
290private:
291 void compute() override { compute<R>();}
292
293 template<class U>
294 typename std::enable_if<!std::is_same<U, TransitionFunction>::value, void>::type
296 {
297 using namespace numeric;
298 auto& result = this->accessValueMutable ();
299 const auto& x0 = accessValueConstCast<T0>(*this->dependency (0));
300 const auto& x1 = accessValueConstCast<T1>(*this->dependency (1));
301 result = cwise(x0) + cwise(x1);
302
303#ifdef DEBUG
304 std::cerr << "=== Add === " << this << std::endl;
305 std::cerr << "x0= " << x0 << std::endl;
306 std::cerr << "x1= " << x1 << std::endl;
307 std::cerr << "res=" << result << std::endl;
308 std::cerr << "=== end Add === " << this << std::endl << std::endl;
309#endif
310 }
311
312 template<class U>
313 typename std::enable_if<std::is_same<U, TransitionFunction>::value, void>::type
315 {
316 using namespace numeric;
317 auto& result = this->accessValueMutable ();
318 const auto& x0 = accessValueConstCast<T0>(*this->dependency (0));
319 const auto& x1 = accessValueConstCast<T1>(*this->dependency (1));
320
321 result = [x0, x1](const VectorLik& x)->VectorLik {return x0(x) + x1(x);};
322 }
323
325};
326
327
328/*************************************************************************
329 * @brief r = x0 - x1 for each component.
330 * - r: R.
331 * - x0: T0.
332 * - x1: T1.
333 *
334 * Values converted to R with the semantics of numeric::convert.
335 * Node construction should be done with the create static method.
336 *
337 * Only defined for N = 2 for now.
338 * The generic version is horrible in C++11 (lack of auto return).
339 * Generic simplification routine is horrible too.
340 */
341template<typename R, typename T0, typename T1> class CWiseSub<R, std::tuple<T0, T1>> : public Value<R>
342{
343public:
344 using Self = CWiseSub;
345
347 static ValueRef<R> create (Context& c, NodeRefVec&& deps, const Dimension<R>& dim)
348 {
349 // Check dependencies
350 checkDependenciesNotNull (typeid (Self), deps);
351 checkDependencyVectorSize (typeid (Self), deps, 2);
352 checkNthDependencyIsValue<T0>(typeid (Self), deps, 0);
353 checkNthDependencyIsValue<T1>(typeid (Self), deps, 1);
354 // Select node implementation
355 bool zeroDep0 = deps[0]->hasNumericalProperty (NumericalProperty::ConstantZero);
356 bool zeroDep1 = deps[1]->hasNumericalProperty (NumericalProperty::ConstantZero);
357 if (zeroDep0 && zeroDep1)
358 {
359 return ConstantZero<R>::create (c, dim);
360 }
361 else if (zeroDep0 && !zeroDep1)
362 {
363 return Convert<R, T1>::create (c, {deps[1]}, dim);
364 }
365 else if (!zeroDep0 && zeroDep1)
366 {
367 return Convert<R, T0>::create (c, {deps[0]}, dim);
368 }
369 else
370 {
371 return cachedAs<Value<R>>(c, std::make_shared<Self>(std::move (deps), dim));
372 }
373 }
374
375 CWiseSub (NodeRefVec&& deps, const Dimension<R>& dim)
376 : Value<R>(std::move (deps)), targetDimension_ (dim) {}
377
378 std::string debugInfo () const override
379 {
380 using namespace numeric;
381 return debug (this->accessValueConst ()) + " targetDim=" + to_string (targetDimension_);
382 }
383
384 // Convert<T> additional arguments = ().
385 bool compareAdditionalArguments (const Node_DF& other) const final
386 {
387 return dynamic_cast<const Self*>(&other) != nullptr;
388 }
389
390 NodeRef derive (Context& c, const Node_DF& node) final
391 {
392 if (&node == this)
393 {
394 return ConstantOne<R>::create (c, targetDimension_);
395 }
396 constexpr std::size_t n = 2;
397 NodeRefVec derivedDeps (n);
398 for (std::size_t i = 0; i < n; ++i)
399 {
400 derivedDeps[i] = this->dependency (i)->derive (c, node);
401 }
402 return Self::create (c, std::move (derivedDeps), targetDimension_);
403 }
404
406 {
407 return Self::create (c, std::move (deps), targetDimension_);
408 }
409
410private:
411 void compute () final
412 {
413 using namespace numeric;
414 auto& result = this->accessValueMutable ();
415 const auto& x0 = accessValueConstCast<T0>(*this->dependency (0));
416 const auto& x1 = accessValueConstCast<T1>(*this->dependency (1));
417 result = cwise(x0) - cwise(x1);
418#ifdef DEBUG
419 std::cerr << "=== Sub === " << this << std::endl;
420 std::cerr << "x0= " << x0 << std::endl;
421 std::cerr << "x1= " << x1 << std::endl;
422 std::cerr << "result= " << result << std::endl;
423 std::cerr << "=== end Sub === " << this << std::endl << std::endl;
424#endif
425 }
426
428};
429
430
431/*************************************************************************
432 * @brief r = sum (x_i), for each component.
433 * - r: R.
434 * - x_i: T.
435 *
436 * Sum of any number of T values into R.
437 * Values converted to R with the semantics of numeric::convert.
438 * Node construction should be done with the create static method.
439 */
440template<typename R, typename T> class CWiseAdd<R, ReductionOf<T>> : public Value<R>
441{
442public:
443 using Self = CWiseAdd;
444
446 static ValueRef<R> create (Context& c, NodeRefVec&& deps, const Dimension<R>& dim)
447 {
448 // Check dependencies
449 checkDependenciesNotNull (typeid (Self), deps);
450 checkDependencyRangeIsValue<T>(typeid (Self), deps, 0, deps.size ());
451 // Remove 0s from deps
452 removeDependenciesIf (deps, [](const NodeRef& dep) -> bool {
453 return dep->hasNumericalProperty (NumericalProperty::ConstantZero);
454 });
455 // Select node implementation
456 if (deps.size () == 0)
457 {
458 return ConstantZero<R>::create (c, dim);
459 }
460 else if (deps.size () == 1)
461 {
462 return Convert<R, T>::create (c, std::move (deps), dim);
463 }
464 else if (deps.size () == 2)
465 {
466 return CWiseAdd<R, std::tuple<T, T>>::create (c, std::move (deps), dim);
467 }
468 else
469 {
470 return cachedAs<Value<R>>(c, std::make_shared<Self>(std::move (deps), dim));
471 }
472 }
473
474 CWiseAdd (NodeRefVec&& deps, const Dimension<R>& dim)
475 : Value<R>(std::move (deps)), targetDimension_ (dim) {}
476
477 std::string debugInfo () const override
478 {
479 using namespace numeric;
480 return debug (this->accessValueConst ()) + " targetDim=" + to_string (targetDimension_);
481 }
482
483 // CWiseAdd additional arguments = ().
484 bool compareAdditionalArguments (const Node_DF& other) const final
485 {
486 return dynamic_cast<const Self*>(&other) != nullptr;
487 }
488
489 NodeRef derive (Context& c, const Node_DF& node) final
490 {
491 if (&node == this)
492 {
494 }
495 const auto n = this->nbDependencies ();
496 NodeRefVec derivedDeps (n);
497 for (std::size_t i = 0; i < n; ++i)
498 {
499 derivedDeps[i] = this->dependency (i)->derive (c, node);
500 }
501 return Self::create (c, std::move (derivedDeps), targetDimension_);
502 }
503
505 {
506 return Self::create (c, std::move (deps), targetDimension_);
507 }
508
509private:
510 void compute() override { compute<T>();}
511
512 template<class U>
513 typename std::enable_if<!std::is_same<U, TransitionFunction>::value, void>::type
515 {
516 using namespace numeric;
517 auto& result = this->accessValueMutable ();
518 result = zero (targetDimension_);
519 for (const auto& depNodeRef : this->dependencies ())
520 {
521 result += accessValueConstCast<T>(*depNodeRef);
522 }
523 }
524
525 template<class U>
526 typename std::enable_if<std::is_same<U, TransitionFunction>::value, void>::type
528 {
529 using namespace numeric;
530 auto& result = this->accessValueMutable ();
531 std::list<const T*> lT;
532 for (const auto& depNodeRef : this->dependencies ())
533 {
534 lT.push_back(&accessValueConstCast<T>(*depNodeRef));
535 }
536
537 result = [lT, this](const VectorLik& x)->VectorLik {
538 VectorLik r = zero (Dimension<VectorLik>(targetDimension_.cols, (Eigen::Index)1));
539
540 for (const auto f:lT)
541 {
542 cwise(r) += cwise((*f)(x));
543 }
544 return r;
545 };
546 }
547
549};
550
551
552/*************************************************************************
553 * @brief r = sum (x_i), for each component either column-wise or row-wise, depending of entry & return types.
554 *
555 * - r: R.
556 * - (x_i): T.
557 *
558 * Sum of any number of values of T into bits of R
559 *.
560 * Node construction should be done with the create static method.
561 */
562
563template<typename R, typename T> class CWiseAdd : public Value<R>
564{
565public:
566 using Self = CWiseAdd;
567
569 static ValueRef<R> create (Context& c, NodeRefVec&& deps, const Dimension<R>& dim)
570 {
571 // Check dependencies
572 checkDependenciesNotNull (typeid (Self), deps);
573 checkDependencyVectorSize (typeid (Self), deps, 1);
574
576 return ConstantZero<R>::create (c, dim);
577 else
578 {
579 return cachedAs<Value<R>>(c, std::make_shared<Self>(std::move (deps), dim));
580 }
581 }
582
583 CWiseAdd (NodeRefVec&& deps, const Dimension<R>& dim)
584 : Value<R>(std::move (deps)), targetDimension_ (dim) {}
585
586 std::string debugInfo () const override
587 {
588 using namespace numeric;
589 return debug (this->accessValueConst ()) + " targetDim=" + to_string (targetDimension_);
590 }
591
592 // CWiseAdd additional arguments = ().
593 bool compareAdditionalArguments (const Node_DF& other) const final
594 {
595 return dynamic_cast<const Self*>(&other) != nullptr;
596 }
597
598 NodeRef derive (Context& c, const Node_DF& node) final
599 {
600 if (&node == this)
601 {
603 }
604 return Self::create (c, {this->dependency(0)->derive (c, node)}, targetDimension_);
605 }
606
608 {
609 return Self::create (c, std::move (deps), targetDimension_);
610 }
611
612private:
613 void compute() override { compute<R, T>();}
614
615 template<class U, class V>
616 typename std::enable_if<(std::is_base_of<V, MatrixLik>::value) && (std::is_same<U, RowLik>::value || std::is_same<U, Eigen::RowVectorXd>::value), void>::type
618 {
619 const auto& mat = accessValueConstCast<T>(*this->dependency(0));
620 this->accessValueMutable () = mat.colwise().sum();
621 }
622
623 template<class U, class V>
624 typename std::enable_if<(std::is_base_of<V, MatrixLik>::value) && (std::is_same<U, VectorLik>::value || std::is_same<U, Eigen::VectorXd>::value), void>::type
626 {
627 const auto& mat = accessValueConstCast<T>(*this->dependency(0));
628 this->accessValueMutable () = mat.rowwise().sum();
629 }
630
631 template<class U, class V>
632 typename std::enable_if<(std::is_base_of<V, VectorLik>::value) || (std::is_same<V, RowLik>::value || std::is_same<V, Eigen::RowVectorXd>::value), void>::type
634 {
635 const auto& vec = accessValueConstCast<T>(*this->dependency(0));
636 this->accessValueMutable () = vec.sum();
637 }
638
639
641};
642
643
644/*************************************************************************
645 * @brief r = sum (p_i * x_i), for each component.
646 * - r: R.
647 * - x_i: T.
648 * - p_i: P
649 *
650 * Sum of any number of T values multiplied per P values into R.
651 * Values converted to R with the semantics of numeric::convert.
652 * Node construction should be done with the create static method.
653 */
654
655template<typename R, typename T, typename P> class CWiseMean<R, ReductionOf<T>, ReductionOf<P>> : public Value<R>
656{
657public:
659
661 // Last dependency is for P component
662 static ValueRef<R> create (Context& c, NodeRefVec&& deps, const Dimension<R>& dim)
663 {
664 // Check dependencies
665 checkDependenciesNotNull (typeid (Self), deps);
666 if (deps.size() % 2 == 1)
667 Exception ("CWiseMean needs an even number of dependencies, got " + std::to_string(deps.size()));
668
669 size_t half = deps.size () / 2;
670
671 checkDependencyRangeIsValue<T>(typeid (Self), deps, 0, half);
672 checkDependencyRangeIsValue<P>(typeid (Self), deps, half, deps.size());
673
674 // Remove both p_i and x_i from deps if one is null
675
676 for (size_t i = 0; i < half; i++)
677 {
678 if (deps[i]->hasNumericalProperty (NumericalProperty::ConstantZero) ||
679 deps[i + half]->hasNumericalProperty (NumericalProperty::ConstantZero))
680 {
681 deps[i] = 0;
682 deps[i + half] = 0;
683 }
684 }
685
686 removeDependenciesIf (deps, [](const NodeRef& dep)->bool {
687 return dep == 0;
688 });
689
690 // if p_i * x_i are all Null
691 if (deps.size() == 0)
692 return ConstantZero<R>::create (c, dim);
693
694 // Select node implementation
695 return cachedAs<Value<R>>(c, std::make_shared<Self>(std::move (deps), dim));
696 }
697
698 CWiseMean (NodeRefVec&& deps, const Dimension<R>& dim)
699 : Value<R>(std::move (deps)), targetDimension_ (dim) {}
700
701 std::string debugInfo () const override
702 {
703 using namespace numeric;
704 return debug (this->accessValueConst ()) + " targetDim=" + to_string (targetDimension_);
705 }
706
707 std::string shape() const override
708 {
709 return "trapezium";
710 }
711
712 std::string color() const override
713 {
714 return "#ffd0d0";
715 }
716
717 std::string description() const override
718 {
719 return "Mean";
720 }
721
722 // CWiseMean additional arguments = ().
723 bool compareAdditionalArguments (const Node_DF& other) const final
724 {
725 return dynamic_cast<const Self*>(&other) != nullptr;
726 }
727
728 NodeRef derive (Context& c, const Node_DF& node) final
729 {
730 if (&node == this)
731 {
732 return ConstantOne<R>::create (c, targetDimension_);
733 }
734 const auto n = this->nbDependencies ();
735 size_t half = n / 2;
736
737 NodeRefVec derivedDeps_T (n);
738 for (std::size_t i = 0; i < half; ++i)
739 {
740 derivedDeps_T[i] = this->dependency (i)->derive (c, node);
741 }
742 for (std::size_t i = half; i < n; ++i)
743 {
744 derivedDeps_T[i] = this->dependency (i);
745 }
746 NodeRef dR_dT = Self::create (c, std::move (derivedDeps_T), targetDimension_);
747
748 NodeRefVec derivedDeps_P(n);
749 for (std::size_t i = 0; i < half; ++i)
750 {
751 derivedDeps_P[i] = this->dependency (i);
752 }
753 for (std::size_t i = half; i < n; ++i)
754 {
755 derivedDeps_P[i] = this->dependency (i)->derive (c, node);
756 }
757
758 NodeRef dR_dP = Self::create (c, std::move (derivedDeps_P), targetDimension_);
759
760 return CWiseAdd<R, std::tuple<R, R>>::create(c, {dR_dT, dR_dP}, targetDimension_);
761 }
762
764 {
765 return Self::create (c, std::move (deps), targetDimension_);
766 }
767
768private:
769 void compute () final
770 {
771 using namespace numeric;
772 auto& result = this->accessValueMutable ();
773 result = zero (targetDimension_);
774 size_t half = this->nbDependencies() / 2;
775 for (size_t i = 0; i < half; i++)
776 {
777 cwise(result) += cwise (accessValueConstCast<P>(*this->dependency(i + half))) * cwise (accessValueConstCast<T>(*this->dependency(i)));
778 }
779 }
780
782};
783
784template<typename R, typename T, typename P> class CWiseMean<R, ReductionOf<T>, P> : public Value<R>
785{
786public:
788
790 // Last dependency is for P component
791 static ValueRef<R> create (Context& c, NodeRefVec&& deps, const Dimension<R>& dim)
792 {
793 // Check dependencies
794 if (deps.size () <= 1)
795 return ConstantZero<R>::create (c, dim);
796 checkDependencyRangeIsValue<T>(typeid (Self), deps, 0, deps.size () - 1);
797 checkNthDependencyIsValue<P>(typeid (Self), deps, deps.size () - 1);
798 // if p_i are all Null
799 if (deps[deps.size() - 1]->hasNumericalProperty (NumericalProperty::ConstantZero))
800 return ConstantZero<R>::create (c, dim);
801 // if x_i are all Null
802 if (std::all_of (deps.begin (), deps.end () - 1, [](const NodeRef& dep)->bool {
803 return dep->hasNumericalProperty (NumericalProperty::ConstantZero);
804 }))
805 {
806 return ConstantZero<R>::create (c, dim);
807 }
808
809 // Select node implementation
810 return cachedAs<Value<R>>(c, std::make_shared<Self>(std::move (deps), dim));
811 }
812
813 CWiseMean (NodeRefVec&& deps, const Dimension<R>& dim)
814 : Value<R>(std::move (deps)), targetDimension_ (dim) {}
815
816 std::string debugInfo () const override
817 {
818 using namespace numeric;
819 return debug (this->accessValueConst ()) + " targetDim=" + to_string (targetDimension_);
820 }
821
822 std::string shape() const override
823 {
824 return "trapezium";
825 }
826
827 std::string color() const override
828 {
829 return "#ffd0d0";
830 }
831
832 std::string description() const override
833 {
834 return "Mean";
835 }
836
837 // CWiseAdd additional arguments = ().
838 bool compareAdditionalArguments (const Node_DF& other) const final
839 {
840 return dynamic_cast<const Self*>(&other) != nullptr;
841 }
842
843 NodeRef derive (Context& c, const Node_DF& node) final
844 {
845 if (&node == this)
846 {
847 return ConstantOne<R>::create (c, targetDimension_);
848 }
849 const auto n = this->nbDependencies ();
850 NodeRefVec derivedDeps_T (n);
851 for (std::size_t i = 0; i < n - 1; ++i)
852 {
853 derivedDeps_T[i] = this->dependency (i)->derive (c, node);
854 }
855 derivedDeps_T[n - 1] = this->dependency (n - 1);
856 NodeRef dR_dT = Self::create (c, std::move (derivedDeps_T), targetDimension_);
857
858 NodeRefVec derivedDeps_P(n);
859 for (std::size_t i = 0; i < n - 1; ++i)
860 {
861 derivedDeps_P[i] = this->dependency (i);
862 }
863 derivedDeps_P[n - 1] = this->dependency (n - 1)->derive (c, node);
864 NodeRef dR_dP = Self::create (c, std::move (derivedDeps_P), targetDimension_);
865 return CWiseAdd<R, std::tuple<R, R>>::create(c, {dR_dT, dR_dP}, targetDimension_);
866 }
867
869 {
870 return Self::create (c, std::move (deps), targetDimension_);
871 }
872
873private:
874 void compute () final
875 {
876 using namespace numeric;
877 auto& result = this->accessValueMutable ();
878 auto& p = accessValueConstCast<P>(*this->dependency(this->nbDependencies() - 1));
879
880#ifdef DEBUG
881 std::cerr << "=== CWiseMean === " << this << std::endl;
882 std::cerr << this->nbDependencies() << std::endl;
883 std::cerr << this->dependency(this->nbDependencies() - 1) << std::endl;
884 std::cerr << "p= " << p << std::endl;
885 std::cerr << "=== end CWiseMean === " << this << std::endl << std::endl;
886#endif
887
888 result = cwise(p)[0] * cwise (accessValueConstCast<T>(*this->dependency(0)));
889 for (Eigen::Index i = 1; i < Eigen::Index(this->nbDependencies() - 1); i++)
890 {
891 cwise(result) += cwise(p)[i] * cwise (accessValueConstCast<T>(*this->dependency(size_t(i))));
892 }
893 }
894
896};
897
898
899/*************************************************************************
900 * @brief r = x0 * x1 for each component.
901 * - r: R.
902 * - x0: T0.
903 * - x1: T1.
904 *
905 * Values converted to R with the semantics of numeric::convert.
906 * Node construction should be done with the create static method.
907 *
908 * Only defined for N = 2 for now (same constraints as CWiseAdd for genericity).
909 */
910
911template<typename R, typename T0, typename T1> class CWiseMul<R, std::tuple<T0, T1>> : public Value<R>
912{
913public:
914 using Self = CWiseMul;
915
917 static ValueRef<R> create (Context& c, NodeRefVec&& deps, const Dimension<R>& dim)
918 {
919 // Check dependencies
920 checkDependenciesNotNull (typeid (Self), deps);
921 checkDependencyVectorSize (typeid (Self), deps, 2);
922 checkNthDependencyIsValue<T0>(typeid (Self), deps, 0);
923 checkNthDependencyIsValue<T1>(typeid (Self), deps, 1);
924 // Return 0 if any 0.
925 if (std::any_of (deps.begin (), deps.end (), [](const NodeRef& dep) -> bool {
926 return dep->hasNumericalProperty (NumericalProperty::ConstantZero);
927 }))
928 {
929 return ConstantZero<R>::create (c, dim);
930 }
931 // Select node implementation
932 bool oneDep0 = deps[0]->hasNumericalProperty (NumericalProperty::ConstantOne);
933 bool oneDep1 = deps[1]->hasNumericalProperty (NumericalProperty::ConstantOne);
934 if (oneDep0 && oneDep1)
935 {
936 return ConstantOne<R>::create (c, dim);
937 }
938 else if (oneDep0 && !oneDep1)
939 {
940 return Convert<R, T1>::create (c, {deps[1]}, dim);
941 }
942 else if (!oneDep0 && oneDep1)
943 {
944 return Convert<R, T0>::create (c, {deps[0]}, dim);
945 }
946 else
947 {
948 return cachedAs<Value<R>>(c, std::make_shared<Self>(std::move (deps), dim));
949 }
950 }
951
952 CWiseMul (NodeRefVec&& deps, const Dimension<R>& dim)
953 : Value<R>(std::move (deps)), targetDimension_ (dim) {}
954
955 std::string debugInfo () const override
956 {
957 using namespace numeric;
958 return debug (this->accessValueConst ()) + " targetDim=" + to_string (targetDimension_);
959 }
960
961 std::string shape() const override
962 {
963 return "house";
964 }
965
966 std::string color() const override
967 {
968 return "#ff9090";
969 }
970
971 std::string description() const override
972 {
973 return "Mul";
974 }
975
976 // CWiseMul additional arguments = ().
977 bool compareAdditionalArguments (const Node_DF& other) const final
978 {
979 return dynamic_cast<const Self*>(&other) != nullptr;
980 }
981
982 NodeRef derive (Context& c, const Node_DF& node) final
983 {
984 if (&node == this)
985 {
986 return ConstantOne<R>::create (c, targetDimension_);
987 }
988 constexpr std::size_t n = 2;
989 NodeRefVec addDeps (n);
990 for (std::size_t i = 0; i < n; ++i)
991 {
992 NodeRefVec ithMulDeps = this->dependencies ();
993 ithMulDeps[i] = this->dependency (i)->derive (c, node);
994 addDeps[i] = Self::create (c, std::move (ithMulDeps), targetDimension_);
995 }
996 return CWiseAdd<R, std::tuple<R, R>>::create (c, std::move (addDeps), targetDimension_);
997 }
998
1000 {
1001 return Self::create (c, std::move (deps), targetDimension_);
1002 }
1003
1004private:
1005 void compute() override { compute<T0, T1>();}
1006
1007 template<class U, class V>
1008 typename std::enable_if<!std::is_same<U, TransitionFunction>::value && !std::is_same<V, TransitionFunction>::value, void>::type
1010 {
1011 using namespace numeric;
1012 auto& result = this->accessValueMutable ();
1013 const auto& x0 = accessValueConstCast<U>(*this->dependency (0));
1014 const auto& x1 = accessValueConstCast<V>(*this->dependency (1));
1015#ifdef DEBUG
1016 std::cerr << "=== Mul === " << this << std::endl;
1017 std::cerr << "x0= " << x0 << std::endl;
1018 std::cerr << "x1= " << x1 << std::endl;
1019#endif
1020 result = cwise (x0) * cwise (x1);
1021
1022#ifdef DEBUG
1023 std::cerr << "result= " << result << std::endl;
1024 std::cerr << "=== end Mul === " << this << std::endl << std::endl;
1025#endif
1026 }
1027
1028 template<class U, class V>
1029 typename std::enable_if<std::is_same<U, TransitionFunction>::value && std::is_same<V, TransitionFunction>::value, void>::type
1031 {
1032 using namespace numeric;
1033 auto& result = this->accessValueMutable ();
1034 const auto& x0 = accessValueConstCast<U>(*this->dependency (0));
1035 const auto& x1 = accessValueConstCast<V>(*this->dependency (1));
1036
1037 result = [x0, x1](const VectorLik& x)->VectorLik {return cwise(x0(x)) * cwise(x1(x));};
1038
1039#ifdef DEBUG
1040 std::cerr << "=== Mul Transition Function X Transition Function === " << this << std::endl;
1041 std::cerr << "=== end Mul Transition Function X Transition Function === " << this << std::endl << std::endl;
1042#endif
1043 }
1044
1045 template<class U, class V>
1046 typename std::enable_if<std::is_same<U, TransitionFunction>::value && !std::is_same<V, TransitionFunction>::value, void>::type
1048 {
1049 using namespace numeric;
1050 auto& result = this->accessValueMutable ();
1051 const auto& x0 = accessValueConstCast<U>(*this->dependency (0));
1052 const auto& x1 = accessValueConstCast<V>(*this->dependency (1));
1053
1054 result = [x0, x1](const VectorLik& x)->VectorLik {return cwise(x0(x)) * cwise(x1);};
1055
1056#ifdef DEBUG
1057 std::cerr << "=== Mul Transition Function X Normal === " << this << std::endl;
1058 std::cerr << "x1= " << x1 << std::endl;
1059 std::cerr << "=== end Mul Transition Function X Normal === " << this << std::endl << std::endl;
1060#endif
1061 }
1062
1063 template<class U, class V>
1064 typename std::enable_if<!std::is_same<U, TransitionFunction>::value && std::is_same<V, TransitionFunction>::value, void>::type
1066 {
1067 using namespace numeric;
1068 auto& result = this->accessValueMutable ();
1069 const auto& x0 = accessValueConstCast<U>(*this->dependency (0));
1070 const auto& x1 = accessValueConstCast<V>(*this->dependency (1));
1071
1072 result = [x0, x1](const VectorLik& x)->VectorLik {return cwise(x1(x)) * cwise(x0);};
1073
1074#ifdef DEBUG
1075 std::cerr << "=== Mul Normal X Transition Function === " << this << std::endl;
1076 std::cerr << "x0= " << x0 << std::endl;
1077 std::cerr << "=== end Mul Normal X Transition Function === " << this << std::endl << std::endl;
1078#endif
1079 }
1080
1082};
1083
1084/*************************************************************************
1085 * @brief r = prod (x_i), for each component.
1086 * - r: R.
1087 * - x_i: T.
1088 *
1089 * Product of any number of T values into R.
1090 * Values converted to R with the semantics of numeric::convert.
1091 * Node construction should be done with the create static method.
1092 */
1093template<typename R, typename T> class CWiseMul<R, ReductionOf<T>> : public Value<R>
1094{
1095public:
1097
1099 static ValueRef<R> create (Context& c, NodeRefVec&& deps, const Dimension<R>& dim)
1100 {
1101 // Check dependencies
1102 checkDependenciesNotNull (typeid (Self), deps);
1103 checkDependencyRangeIsValue<T>(typeid (Self), deps, 0, deps.size ());
1104 // If there is a 0 return 0.
1105 if (std::any_of (deps.begin (), deps.end (), [](const NodeRef& dep) -> bool {
1106 return dep->hasNumericalProperty (NumericalProperty::ConstantZero);
1107 }))
1108 {
1109 return ConstantZero<R>::create (c, dim);
1110 }
1111 // Remove 1s from deps
1112 removeDependenciesIf (deps, [](const NodeRef& dep)->bool {
1113 return dep->hasNumericalProperty (NumericalProperty::ConstantOne);
1114 });
1115 // Select node implementation
1116 if (deps.size () == 0)
1117 {
1118 return ConstantOne<R>::create (c, dim);
1119 }
1120 else if (deps.size () == 1)
1121 {
1122 return Convert<R, T>::create (c, std::move (deps), dim);
1123 }
1124 else if (deps.size () == 2)
1125 {
1126 return CWiseMul<R, std::tuple<T, T>>::create (c, std::move (deps), dim);
1127 }
1128 else
1129 {
1130 return cachedAs<Value<R>>(c, std::make_shared<Self>(std::move (deps), dim));
1131 }
1132 }
1133
1134 CWiseMul (NodeRefVec&& deps, const Dimension<R>& dim)
1135 : Value<R>(std::move (deps)), targetDimension_ (dim) {}
1136
1137 std::string debugInfo () const override
1138 {
1139 using namespace numeric;
1140 return debug (this->accessValueConst ()) + " targetDim=" + to_string (targetDimension_);
1141 }
1142
1143 // CWiseMul additional arguments = ().
1144 bool compareAdditionalArguments (const Node_DF& other) const final
1145 {
1146 return dynamic_cast<const Self*>(&other) != nullptr;
1147 }
1148
1149 NodeRef derive (Context& c, const Node_DF& node) final
1150 {
1151 if (&node == this)
1152 {
1153 return ConstantOne<R>::create (c, targetDimension_);
1154 }
1155 const auto n = this->nbDependencies ();
1156 NodeRefVec addDeps (n);
1157 for (std::size_t i = 0; i < n; ++i)
1158 {
1159 NodeRefVec ithMulDeps = this->dependencies ();
1160 ithMulDeps[i] = this->dependency (i)->derive (c, node);
1161 addDeps[i] = Self::create (c, std::move (ithMulDeps), targetDimension_);
1162 }
1163 return CWiseAdd<R, ReductionOf<R>>::create (c, std::move (addDeps), targetDimension_);
1164 }
1165
1167 {
1168 return Self::create (c, std::move (deps), targetDimension_);
1169 }
1170
1171private:
1172 void compute () final
1173 {
1174 using namespace numeric;
1175 auto& result = this->accessValueMutable ();
1176 result = one (targetDimension_);
1177 for (const auto& depNodeRef : this->dependencies ())
1178 {
1179 cwise(result) *= cwise (accessValueConstCast<T>(*depNodeRef));
1180 }
1181 }
1182
1184};
1185
1186
1187/*************************************************************************
1188 * @brief r = x0 / x1 for each component.
1189 * - r: R.
1190 * - x0: T0.
1191 * - x1: T1.
1192 *
1193 * Values converted to R with the semantics of numeric::convert.
1194 * Node construction should be done with the create static method.
1195 *
1196 * Only defined for N = 2 for now (same constraints as CWiseAdd for genericity).
1197 */
1198template<typename R, typename T0, typename T1> class CWiseDiv<R, std::tuple<T0, T1>> : public Value<R>
1199{
1200public:
1202
1204 static ValueRef<R> create (Context& c, NodeRefVec&& deps, const Dimension<R>& dim)
1205 {
1206 // Check dependencies
1207 checkDependenciesNotNull (typeid (Self), deps);
1208 checkDependencyVectorSize (typeid (Self), deps, 2);
1209 checkNthDependencyIsValue<T0>(typeid (Self), deps, 0);
1210 checkNthDependencyIsValue<T1>(typeid (Self), deps, 1);
1211 // Select node
1212 if (deps[1]->hasNumericalProperty (NumericalProperty::ConstantOne))
1213 return Convert<R, T0>::create (c, {deps[0]}, dim);
1214 else if (deps[0]->hasNumericalProperty (NumericalProperty::ConstantOne))
1215 return CWiseInverse<R>::create (c, {deps[1]}, dim);
1216 else
1217 return cachedAs<Value<R>>(c, std::make_shared<Self>(std::move (deps), dim));
1218 }
1219
1220 CWiseDiv (NodeRefVec&& deps, const Dimension<R>& dim)
1221 : Value<R>(std::move (deps)), targetDimension_ (dim) {}
1222
1223 std::string debugInfo () const override
1224 {
1225 using namespace numeric;
1226 return debug (this->accessValueConst ()) + " targetDim=" + to_string (targetDimension_);
1227 }
1228
1229 std::string shape() const override
1230 {
1231 return "house";
1232 }
1233
1234 std::string color() const override
1235 {
1236 return "#ff9090";
1237 }
1238
1239 std::string description() const override
1240 {
1241 return "Div";
1242 }
1243
1244 // CWiseDiv additional arguments = ().
1245 bool compareAdditionalArguments (const Node_DF& other) const final
1246 {
1247 return dynamic_cast<const Self*>(&other) != nullptr;
1248 }
1249
1250 NodeRef derive (Context& c, const Node_DF& node) final
1251 {
1252 if (&node == this)
1253 {
1254 return ConstantOne<R>::create (c, targetDimension_);
1255 }
1256
1257 auto f0 = this->dependency (0);
1258 auto f1 = this->dependency (1);
1259 auto fp0 = this->dependency (0)->derive (c, node);
1260 auto fp1 = this->dependency (1)->derive (c, node);
1261
1262 auto upv = CWiseMul<R, std::tuple<T0, T1>>::create(c, {fp0, f1}, targetDimension_);
1263 auto vpu = CWiseMul<R, std::tuple<T0, T1>>::create(c, {fp1, f0}, targetDimension_);
1264 auto diff = CWiseSub<R, std::tuple<R, R>>::create(c, {upv, vpu}, targetDimension_);
1265
1266 return CWiseMul<R, std::tuple<R, R>>::create (
1267 c, {CWiseConstantPow<R>::create(c, {f1}, -2., 1., targetDimension_), diff}, targetDimension_);
1268 }
1269
1271 {
1272 return Self::create (c, std::move (deps), targetDimension_);
1273 }
1274
1275private:
1276 void compute() override { compute<T0, T1>();}
1277
1278 template<class U, class V>
1279 typename std::enable_if<!std::is_same<U, TransitionFunction>::value && !std::is_same<V, TransitionFunction>::value, void>::type
1281 {
1282 using namespace numeric;
1283 auto& result = this->accessValueMutable ();
1284 const auto& x0 = accessValueConstCast<U>(*this->dependency (0));
1285 const auto& x1 = accessValueConstCast<V>(*this->dependency (1));
1286 result = cwise (x0) / cwise (x1);
1287 }
1288
1289 template<class U, class V>
1290 typename std::enable_if<std::is_same<U, TransitionFunction>::value && std::is_same<V, TransitionFunction>::value, void>::type
1292 {
1293 using namespace numeric;
1294 auto& result = this->accessValueMutable ();
1295 const auto& x0 = accessValueConstCast<U>(*this->dependency (0));
1296 const auto& x1 = accessValueConstCast<V>(*this->dependency (1));
1297
1298 result = [x0, x1](const VectorLik& x)->VectorLik {return cwise(x0(x)) / cwise(x1(x));};
1299 }
1300
1301 template<class U, class V>
1302 typename std::enable_if<std::is_same<U, TransitionFunction>::value && !std::is_same<V, TransitionFunction>::value, void>::type
1304 {
1305 using namespace numeric;
1306 auto& result = this->accessValueMutable ();
1307 const auto& x0 = accessValueConstCast<U>(*this->dependency (0));
1308 const auto& x1 = accessValueConstCast<V>(*this->dependency (1));
1309
1310 result = [x0, x1](const VectorLik& x)->VectorLik {return cwise(x0(x)) / cwise(x1);};
1311 }
1312
1313 template<class U, class V>
1314 typename std::enable_if<!std::is_same<U, TransitionFunction>::value && std::is_same<V, TransitionFunction>::value, void>::type
1316 {
1317 using namespace numeric;
1318 auto& result = this->accessValueMutable ();
1319 const auto& x0 = accessValueConstCast<U>(*this->dependency (0));
1320 const auto& x1 = accessValueConstCast<V>(*this->dependency (1));
1321
1322 result = [x0, x1](const VectorLik& x)->VectorLik {return cwise(x1(x)) / cwise(x0);};
1323 }
1324
1326};
1327
1328
1329/*************************************************************************
1330 * @brief r = -x, for each component.
1331 * - r, x: T.
1332 *
1333 * Node construction should be done with the create static method.
1334 */
1335
1336template<typename T> class CWiseNegate : public Value<T>
1337{
1338public:
1340
1342 static ValueRef<T> create (Context& c, NodeRefVec&& deps, const Dimension<T>& dim)
1343 {
1344 // Check dependencies
1345 checkDependenciesNotNull (typeid (Self), deps);
1346 checkDependencyVectorSize (typeid (Self), deps, 1);
1347 checkNthDependencyIsValue<T>(typeid (Self), deps, 0);
1348 // Select node
1350 {
1351 return ConstantZero<T>::create (c, dim);
1352 }
1353 else
1354 {
1355 return cachedAs<Value<T>>(c, std::make_shared<Self>(std::move (deps), dim));
1356 }
1357 }
1358
1360 : Value<T>(std::move (deps)), targetDimension_ (dim) {}
1361
1362 std::string debugInfo () const override
1363 {
1364 using namespace numeric;
1365 return debug (this->accessValueConst ()) + " targetDim=" + to_string (targetDimension_);
1366 }
1367
1368 // CWiseNegate additional arguments = ().
1369 bool compareAdditionalArguments (const Node_DF& other) const final
1370 {
1371 return dynamic_cast<const Self*>(&other) != nullptr;
1372 }
1373
1374 NodeRef derive (Context& c, const Node_DF& node) final
1375 {
1376 if (&node == this)
1377 {
1379 }
1380 return Self::create (c, {this->dependency (0)->derive (c, node)}, targetDimension_);
1381 }
1382
1384 {
1385 return Self::create (c, std::move (deps), targetDimension_);
1386 }
1387
1388private:
1389 void compute () final
1390 {
1391 using namespace numeric;
1392 auto& result = this->accessValueMutable ();
1393 const auto& x = accessValueConstCast<T>(*this->dependency (0));
1394 result = -x;
1395 }
1396
1398};
1399
1400
1401/*************************************************************************
1402 * @brief r = 1/x for each component.
1403 * - r, x: T.
1404 *
1405 * Node construction should be done with the create static method.
1406 */
1407
1408template<typename T> class CWiseInverse : public Value<T>
1409{
1410public:
1412
1414 static ValueRef<T> create (Context& c, NodeRefVec&& deps, const Dimension<T>& dim)
1415 {
1416 // Check dependencies
1417 checkDependenciesNotNull (typeid (Self), deps);
1418 checkDependencyVectorSize (typeid (Self), deps, 1);
1419 checkNthDependencyIsValue<T>(typeid (Self), deps, 0);
1420 // Select node
1422 {
1423 return ConstantOne<T>::create (c, dim);
1424 }
1425 else
1426 {
1427 return cachedAs<Value<T>>(c, std::make_shared<Self>(std::move (deps), dim));
1428 }
1429 }
1430
1432 : Value<T>(std::move (deps)), targetDimension_ (dim) {}
1433
1434 std::string debugInfo () const override
1435 {
1436 using namespace numeric;
1437 return debug (this->accessValueConst ()) + " targetDim=" + to_string (targetDimension_);
1438 }
1439
1440 // CWiseInverse additional arguments = ().
1441 bool compareAdditionalArguments (const Node_DF& other) const final
1442 {
1443 return dynamic_cast<const Self*>(&other) != nullptr;
1444 }
1445
1446 NodeRef derive (Context& c, const Node_DF& node) final
1447 {
1448 if (&node == this)
1449 {
1451 }
1452 // -1/x^2 * x'
1453 const auto& dep = this->dependency (0);
1455 c, {CWiseConstantPow<T>::create (c, {dep}, -2., -1., targetDimension_), dep->derive (c, node)},
1457 }
1458
1460 {
1461 return Self::create (c, std::move (deps), targetDimension_);
1462 }
1463
1464private:
1465 void compute () final
1466 {
1467 using namespace numeric;
1468 auto& result = this->accessValueMutable ();
1469 const auto& x = accessValueConstCast<T>(*this->dependency (0));
1470 result = inverse (cwise (x));
1471 }
1472
1474};
1475
1476
1477/*************************************************************************
1478 * @brief r = log(x) for each component.
1479 * - r, x: T.
1480 *
1481 * Node construction should be done with the create static method.
1482 */
1483
1484using std::log;
1485
1486template<typename T> class CWiseLog : public Value<T>
1487{
1488public:
1490
1492 static ValueRef<T> create (Context& c, NodeRefVec&& deps, const Dimension<T>& dim)
1493 {
1494 // Check dependencies
1495 checkDependenciesNotNull (typeid (Self), deps);
1496 checkDependencyVectorSize (typeid (Self), deps, 1);
1497 checkNthDependencyIsValue<T>(typeid (Self), deps, 0);
1498 // Select node
1500 {
1501 return ConstantZero<T>::create (c, dim);
1502 }
1503 else
1504 {
1505 return cachedAs<Value<T>>(c, std::make_shared<Self>(std::move (deps), dim));
1506 }
1507 }
1508
1509 CWiseLog (NodeRefVec&& deps, const Dimension<T>& dim)
1510 : Value<T>(std::move (deps)), targetDimension_ (dim) {}
1511
1512 std::string debugInfo () const override
1513 {
1514 using namespace numeric;
1515 return debug (this->accessValueConst ()) + " targetDim=" + to_string (targetDimension_);
1516 }
1517
1518 // CWiseLog additional arguments = ().
1519 bool compareAdditionalArguments (const Node_DF& other) const final
1520 {
1521 return dynamic_cast<const Self*>(&other) != nullptr;
1522 }
1523
1524 NodeRef derive (Context& c, const Node_DF& node) final
1525 {
1526 if (&node == this)
1527 {
1529 }
1530 // x'/x
1531 const auto& dep = this->dependency (0);
1533 c, {CWiseInverse<T>::create (c, {dep}, targetDimension_), dep->derive (c, node)}, targetDimension_);
1534 }
1535
1537 {
1538 return Self::create (c, std::move (deps), targetDimension_);
1539 }
1540
1541private:
1542 void compute () final
1543 {
1544 using namespace numeric;
1545 auto& result = this->accessValueMutable ();
1546 const auto& x = accessValueConstCast<T>(*this->dependency (0));
1547 result = log (cwise (x));
1548 }
1549
1551};
1552
1553
1554/*************************************************************************
1555 * @brief r = exp(x) for each component.
1556 * - r, x: T.
1557 *
1558 * Node construction should be done with the create static method.
1559 */
1560
1561using std::exp;
1562
1563template<typename T> class CWiseExp : public Value<T>
1564{
1565public:
1567
1569 static ValueRef<T> create (Context& c, NodeRefVec&& deps, const Dimension<T>& dim)
1570 {
1571 // Check dependencies
1572 checkDependenciesNotNull (typeid (Self), deps);
1573 checkDependencyVectorSize (typeid (Self), deps, 1);
1574 checkNthDependencyIsValue<T>(typeid (Self), deps, 0);
1575 // Select node
1577 {
1578 return ConstantOne<T>::create (c, dim);
1579 }
1580 else
1581 {
1582 return cachedAs<Value<T>>(c, std::make_shared<Self>(std::move (deps), dim));
1583 }
1584 }
1585
1586 CWiseExp (NodeRefVec&& deps, const Dimension<T>& dim)
1587 : Value<T>(std::move (deps)), targetDimension_ (dim) {}
1588
1589 std::string debugInfo () const override
1590 {
1591 using namespace numeric;
1592 return debug (this->accessValueConst ()) + " targetDim=" + to_string (targetDimension_);
1593 }
1594
1595 // CWiseExp additional arguments = ().
1596 bool compareAdditionalArguments (const Node_DF& other) const final
1597 {
1598 return dynamic_cast<const Self*>(&other) != nullptr;
1599 }
1600
1601 NodeRef derive (Context& c, const Node_DF& node) final
1602 {
1603 if (&node == this)
1604 {
1606 }
1607 // x'* exp(x)
1608 const auto& dep = this->dependency (0);
1610 c, {this->shared_from_this(), dep->derive (c, node)}, targetDimension_);
1611 }
1612
1614 {
1615 return Self::create (c, std::move (deps), targetDimension_);
1616 }
1617
1618private:
1619 void compute () final
1620 {
1621 using namespace numeric;
1622 auto& result = this->accessValueMutable ();
1623 const auto& x = accessValueConstCast<T>(*this->dependency (0));
1624 result = exp (cwise (x));
1625 }
1626
1628};
1629
1630
1631/*************************************************************************
1632 * @brief r = factor * pow (x, exponent) for each component.
1633 * - r, x: T.
1634 * - exponent, factor: double (constant parameter of the node).
1635 *
1636 * Node construction should be done with the create static method.
1637 */
1638
1639template<typename T> class CWiseConstantPow : public Value<T>
1640{
1641public:
1643
1645 static ValueRef<T> create (Context& c, NodeRefVec&& deps, double exponent, double factor,
1646 const Dimension<T>& dim)
1647 {
1648 // Check dependencies
1649 checkDependenciesNotNull (typeid (Self), deps);
1650 checkDependencyVectorSize (typeid (Self), deps, 1);
1651 checkNthDependencyIsValue<T>(typeid (Self), deps, 0);
1652 // Select node implementation
1653 if (exponent == 0. || deps[0]->hasNumericalProperty (NumericalProperty::ConstantOne))
1654 {
1655 // pow (x, exponent) == 1
1656 using namespace numeric;
1657 return NumericConstant<T>::create (c, factor * one (dim));
1658 }
1659 else if (exponent == 1.)
1660 {
1661 // pow (x, exponent) == x
1663 c, {NumericConstant<double>::create (c, factor), deps[0]}, dim);
1664 }
1665 else if (exponent == -1.)
1666 {
1667 // pow (x, exponent) = 1/x
1669 c,
1670 {NumericConstant<double>::create (c, factor), CWiseInverse<T>::create (c, std::move (deps), dim)},
1671 dim);
1672 }
1673 else
1674 {
1675 return cachedAs<Value<T>>(c, std::make_shared<Self>(std::move (deps), exponent, factor, dim));
1676 }
1677 }
1678
1679 CWiseConstantPow (NodeRefVec&& deps, double exponent, double factor, const Dimension<T>& dim)
1680 : Value<T>(std::move (deps)), targetDimension_ (dim), exponent_ (exponent), factor_ (factor) {}
1681
1682 std::string debugInfo () const override
1683 {
1684 using namespace numeric;
1685 return debug (this->accessValueConst ()) + " targetDim=" + to_string (targetDimension_) +
1686 " exponent=" + std::to_string (exponent_) + " factor=" + std::to_string (factor_);
1687 }
1688
1689 // CWiseConstantPow additional arguments = (exponent_, factor_).
1690 bool compareAdditionalArguments (const Node_DF& other) const final
1691 {
1692 const auto* derived = dynamic_cast<const Self*>(&other);
1693 return derived != nullptr && exponent_ == derived->exponent_ && factor_ == derived->factor_;
1694 }
1695 std::size_t hashAdditionalArguments () const final
1696 {
1697 std::size_t seed = 0;
1698 combineHash (seed, exponent_);
1699 combineHash (seed, factor_);
1700 return seed;
1701 }
1702
1703 NodeRef derive (Context& c, const Node_DF& node) final
1704 {
1705 if (&node == this)
1706 {
1708 }
1709 // factor * (exponent * x^(exponent - 1)) * x'
1710 const auto& dep = this->dependency (0);
1711 auto dpow = Self::create (c, {dep}, exponent_ - 1., factor_ * exponent_, targetDimension_);
1712 return CWiseMul<T, std::tuple<T, T>>::create (c, {dpow, dep->derive (c, node)}, targetDimension_);
1713 }
1714
1716 {
1717 return Self::create (c, std::move (deps), exponent_, factor_, targetDimension_);
1718 }
1719
1720private:
1721 void compute () final
1722 {
1723 using namespace numeric;
1724 auto& result = this->accessValueMutable ();
1725 const auto& x = accessValueConstCast<T>(*this->dependency (0));
1726 result = factor_ * pow (cwise (x), exponent_);
1727 }
1728
1731 double factor_;
1732};
1733
1734
1735/*************************************************************************
1736 * @brief r = x0 * x1 (dot product).
1737 * - r: double.
1738 * - x0: T0 (vector-like).
1739 * - x1: T1 (vector-like).
1740 *
1741 * Node construction should be done with the create static method.
1742 */
1743template<typename R, typename T0, typename T1> class ScalarProduct : public Value<R>
1744{
1745public:
1747
1750 {
1751 // Check dependencies
1752 checkDependenciesNotNull (typeid (Self), deps);
1753 checkDependencyVectorSize (typeid (Self), deps, 2);
1754 checkNthDependencyIsValue<T0>(typeid (Self), deps, 0);
1755 checkNthDependencyIsValue<T1>(typeid (Self), deps, 1);
1756 // Select node
1759 {
1761 }
1762 else
1763 {
1764 return cachedAs<Value<R>>(c, std::make_shared<Self>(std::move (deps)));
1765 }
1766 }
1767
1768 ScalarProduct (NodeRefVec&& deps) : Value<R>(std::move (deps)) {}
1769
1770 std::string debugInfo () const override
1771 {
1772 using namespace numeric;
1773 return debug (this->accessValueConst ());
1774 }
1775
1776 // ScalarProduct additional arguments = ().
1777 bool compareAdditionalArguments (const Node_DF& other) const final
1778 {
1779 return dynamic_cast<const Self*>(&other) != nullptr;
1780 }
1781
1782 NodeRef derive (Context& c, const Node_DF& node) final
1783 {
1784 if (&node == this)
1785 {
1786 return ConstantOne<R>::create (c, Dimension<R> ());
1787 }
1788 const auto& x0 = this->dependency (0);
1789 const auto& x1 = this->dependency (1);
1790 auto dx0_prod = Self::create (c, {x0->derive (c, node), x1});
1791 auto dx1_prod = Self::create (c, {x0, x1->derive (c, node)});
1792 return CWiseAdd<R, std::tuple<R, R>>::create (c, {dx0_prod, dx1_prod},
1793 Dimension<R> ());
1794 }
1795
1797 {
1798 return Self::create (c, std::move (deps));
1799 }
1800
1801private:
1802 void compute () final
1803 {
1804 auto& result = this->accessValueMutable ();
1805 const auto& x0 = accessValueConstCast<T0>(*this->dependency (0));
1806 const auto& x1 = accessValueConstCast<T1>(*this->dependency (1));
1807 auto d = x0.dot (x1); // Using lhs.dot(rhs) method from Eigen only
1808 result = d;
1809 }
1810};
1811
1812
1813/*************************************************************************
1814 * @brief r = sum_{v in m} log (v).
1815 * - r: double.
1816 * - m: F (matrix-like type).
1817 *
1818 * or r = sum_{i in 0:len(m)-1} log (p[i]*m[i]).
1819 * - r: DataLik
1820 * - m: F (matrix-like type).
1821 * - p: <Eigen::RowVectorXd>
1822 *
1823 * The node has no dimension (double).
1824 * The dimension of m should be provided for derivation.
1825 * Node construction should be done with the create static method.
1826 */
1827
1828template<typename F> class SumOfLogarithms : public Value<DataLik>
1829{
1830public:
1832
1834 static ValueRef<DataLik> create (Context& c, NodeRefVec&& deps, const Dimension<F>& mDim)
1835 {
1836 checkDependenciesNotNull (typeid (Self), deps);
1837 checkDependencyVectorMinSize (typeid (Self), deps, 1);
1838 checkNthDependencyIsValue<F>(typeid (Self), deps, 0);
1839 if (deps.size() == 2)
1840 checkNthDependencyIsValue<Eigen::RowVectorXi>(typeid (Self), deps, 1);
1841 return cachedAs<Value<DataLik>>(c, std::make_shared<Self>(std::move (deps), mDim));
1842 }
1843
1845 : Value<DataLik>(std::move (deps)), mTargetDimension_ (mDim) // , temp_() {
1846 // if (dependencies().size()==2)
1847 // {
1848 // const auto & p = accessValueConstCast<Eigen::VectorXi> (*this->dependency (1));
1849 // temp_.resize(p.size());
1850 // }
1851 {}
1852
1853 std::string debugInfo () const override
1854 {
1855 using namespace numeric;
1856 return debug (this->accessValueConst ());
1857 }
1858
1859 // SumOfLogarithms additional arguments = ().
1860 bool compareAdditionalArguments (const Node_DF& other) const final
1861 {
1862 return dynamic_cast<const Self*>(&other) != nullptr;
1863 }
1864
1865 NodeRef derive (Context& c, const Node_DF& node) final
1866 {
1867 const auto& m = this->dependency (0);
1868 auto dm_dn = m->derive (c, node);
1869 auto m_inverse = CWiseInverse<F>::create (c, {m}, mTargetDimension_);
1870 if (nbDependencies() == 2 && dependency(1))
1871 {
1873 return ScalarProduct<typename F::Scalar, F, F>::create (c, {std::move (dm_dn), std::move (m2)});
1874 }
1875 else
1876 return ScalarProduct<typename F::Scalar, F, F>::create (c, {std::move (dm_dn), std::move (m_inverse)});
1877 }
1878
1880 {
1881 return Self::create (c, std::move (deps), mTargetDimension_);
1882 }
1883
1884private:
1885 void compute() override { compute<F>();}
1886
1887 template<class G = F>
1888 typename std::enable_if<std::is_convertible<G*, Eigen::DenseBase<G>*>::value, void>::type
1890 {
1891#ifdef DEBUG
1892 std::cerr << "=== SumOfLogarithms === " << this << std::endl;
1893#endif
1894
1895 auto& result = this->accessValueMutable ();
1896
1897 const auto& m = accessValueConstCast<F>(*this->dependency (0));
1898
1899 if (nbDependencies() == 1)
1900 {
1901 const ExtendedFloat product = m.unaryExpr ([](double d) {
1902 ExtendedFloat ef{d};
1903 ef.normalize_small ();
1904 return ef;
1905 }).redux ([](const ExtendedFloat& lhs, const ExtendedFloat& rhs) {
1906 auto r = ExtendedFloat::denorm_mul (lhs, rhs);
1907 r.normalize_small ();
1908 return r;
1909 });
1910 result = product.log();
1911#ifdef DEBUG
1912 std::cerr << "product= " << product << std::endl;
1913 std::cerr << "result log= " << result << std::endl;
1914#endif
1915 }
1916 else
1917 {
1918 const auto& p = accessValueConstCast<Eigen::RowVectorXi>(*this->dependency (1));
1919 // Old version:
1920 // double resold = (numeric::cwise(m).log() * numeric::cwise(p)).sum();
1921
1922 temp_ = m.unaryExpr ([](double d) {
1923 ExtendedFloat ef{d};
1924 ef.normalize ();
1925 return ef;
1926 });
1927
1928 for (Eigen::Index i = 0; i < Eigen::Index(p.size()); i++)
1929 {
1930 temp_[i] = temp_[i].pow(p[i]);
1931 }
1932
1933 const ExtendedFloat product = temp_.redux ([](const ExtendedFloat& lhs, const ExtendedFloat& rhs) {
1934 auto r = ExtendedFloat::denorm_mul (lhs, rhs);
1935 r.normalize ();
1936 return r;
1937 });
1938
1939 result = product.log ();
1940#ifdef DEBUG
1941 std::cerr << "PRODUCT= " << product << std::endl;
1942 std::cerr << "RESULT log= " << result << std::endl;
1943#endif
1944 }
1945#ifdef DEBUG
1946 std::cerr << "=== end SumOfLogarithms === " << this << std::endl;
1947#endif
1948 }
1949
1950 template<class G = F>
1951 typename std::enable_if<std::is_convertible<G*, ExtendedFloatEigenBase<G>*>::value, void>::type
1953 {
1954#ifdef DEBUG
1955 std::cerr << "=== SumOfLogarithms === " << this << std::endl;
1956#endif
1957
1958 auto& result = this->accessValueMutable ();
1959
1960 const auto& m = accessValueConstCast<F>(*this->dependency (0));
1961
1962 if (nbDependencies() == 1)
1963 {
1964 const ExtendedFloat product = m.float_part().unaryExpr ([](double d) {
1965 ExtendedFloat ef{d};
1966 ef.normalize ();
1967 return ef;
1968 }).redux ([](const ExtendedFloat& lhs, const ExtendedFloat& rhs) {
1969 auto r = ExtendedFloat::denorm_mul (lhs, rhs);
1970 r.normalize_small ();
1971 return r;
1972 });
1973 result = product.log() + double(m.exponent_part() * m.size()) * ExtendedFloat::ln_radix;
1974#ifdef DEBUG
1975 std::cerr << "product= " << product << "* 2^" << m.size() * m.exponent_part() << std::endl;
1976 std::cerr << "result log= " << result << std::endl;
1977#endif
1978 }
1979 else
1980 {
1981 const auto& p = accessValueConstCast<Eigen::RowVectorXi>(*this->dependency (1));
1982 // Old version:
1983 // double resold = (numeric::cwise(m).log() * numeric::cwise(p)).sum();
1984
1985 temp_ = m.float_part().unaryExpr ([](double d) {
1986 ExtendedFloat ef{d};
1987 ef.normalize ();
1988 return ef;
1989 });
1990
1991 for (Eigen::Index i = 0; i < Eigen::Index(p.size()); i++)
1992 {
1993 temp_[i] = temp_[i].pow(p[i]);
1994 }
1995
1996 const ExtendedFloat product = temp_.redux ([](const ExtendedFloat& lhs, const ExtendedFloat& rhs) {
1997 auto r = ExtendedFloat::denorm_mul (lhs, rhs);
1998 r.normalize ();
1999 return r;
2000 });
2001
2002 result = product.log () + double(m.exponent_part() * p.sum()) * ExtendedFloat::ln_radix;
2003#ifdef DEBUG
2004 std::cerr << "RESULT log= " << result << std::endl;
2005#endif
2006 }
2007#ifdef DEBUG
2008 std::cerr << "=== end SumOfLogarithms === " << this << std::endl;
2009#endif
2010 }
2011
2013
2014 Eigen::Matrix<ExtendedFloat, 1, Eigen::Dynamic> temp_;
2015};
2016
2017
2018/*************************************************************************
2019 * @brief r = log(sum_i p_i * exp (v_i))
2020 * - r: R.
2021 * - v: T0 (vector like type).
2022 * - p: T1 (vector like type).
2023 *
2024 * The node has no dimension (double).
2025 * Node construction should be done with the create static method.
2026 */
2027
2028template<typename R, typename T0, typename T1> class LogSumExp : public Value<R>
2029{
2030public:
2032
2034 static ValueRef<R> create (Context& c, NodeRefVec&& deps, const Dimension<T0>& mDim)
2035 {
2036 checkDependenciesNotNull (typeid (Self), deps);
2037 checkDependencyVectorSize (typeid (Self), deps, 2);
2038 checkNthDependencyIsValue<T0>(typeid (Self), deps, 0);
2039 checkNthDependencyIsValue<T1>(typeid (Self), deps, 1);
2040 // Select node
2041 return cachedAs<Value<R>>(c, std::make_shared<Self>(std::move (deps), mDim));
2042 }
2043
2044 LogSumExp (NodeRefVec&& deps, const Dimension<T0>& mDim)
2045 : Value<R>(std::move (deps)), mTargetDimension_ (mDim) {}
2046
2047 std::string debugInfo () const override
2048 {
2049 using namespace numeric;
2050 return debug (this->accessValueConst ());
2051 }
2052
2053 // LogSumExp additional arguments = ().
2054 bool compareAdditionalArguments (const Node_DF& other) const final
2055 {
2056 return dynamic_cast<const Self*>(&other) != nullptr;
2057 }
2058
2059 NodeRef derive (Context& c, const Node_DF& node) final
2060 {
2061 if (&node == this)
2062 {
2064 }
2065 const auto& v = this->dependency (0);
2066 const auto& p = this->dependency (1);
2067 auto diffvL = CWiseSub<T0, std::tuple<R, T0>>::create(c, {this->shared_from_this(), v}, mTargetDimension_);
2068 auto expdiffvL = CWiseExp<T0>::create(c, {std::move(diffvL)}, mTargetDimension_);
2069 auto dp_prod = ScalarProduct<R, T0, T1>::create (c, {expdiffvL, p->derive (c, node)});
2070 auto pexpdiffvL = CWiseMul<T0, std::tuple<T0, T1>>::create(c, {expdiffvL, p}, mTargetDimension_);
2071 auto dv_prod = ScalarProduct<R, T0, T1>::create (c, {std::move(pexpdiffvL), v->derive (c, node)});
2072 return CWiseAdd<R, std::tuple<R, R>>::create (c, {std::move (dv_prod), std::move (dp_prod)}, Dimension<R>());
2073 }
2074
2076 {
2077 return Self::create (c, std::move (deps), mTargetDimension_);
2078 }
2079
2080private:
2081 void compute () final
2082 {
2083 using namespace numeric;
2084 auto& result = this->accessValueMutable ();
2085 const auto& v = accessValueConstCast<T0>(*this->dependency (0));
2086 const auto& p = accessValueConstCast<T1>(*this->dependency (1));
2087
2088 auto M = v.maxCoeff();
2089 if (isinf(M))
2090 result = M;
2091 else
2092 {
2093 T0 v2;
2094 v2 = exp(cwise(v - T0::Constant(mTargetDimension_.rows, mTargetDimension_.cols, M)));
2095 auto ve = v2.dot(p);
2096 result = log(ve) + M;
2097 }
2098 }
2099
2101};
2102
2103
2104/*************************************************************************
2105 * @brief r = x0 * x1 (matrix product).
2106 * - r: R (matrix).
2107 * - x0: T0 (matrix), allows NumericalDependencyTransform.
2108 * - x1: T1 (matrix), allows NumericalDependencyTransform.
2109 *
2110 * Node construction should be done with the create static method.
2111 */
2112template<typename R, typename T0, typename T1> class MatrixProduct : public Value<R>
2113{
2114public:
2118
2120 static ValueRef<R> create (Context& c, NodeRefVec&& deps, const Dimension<R>& dim)
2121 {
2122 // Check dependencies
2123 checkDependenciesNotNull (typeid (Self), deps);
2124 checkDependencyVectorSize (typeid (Self), deps, 2);
2125 checkNthDependencyIsValue<DepT0>(typeid (Self), deps, 0);
2126 checkNthDependencyIsValue<DepT1>(typeid (Self), deps, 1);
2127 // Return 0 if any 0.
2128 if (std::any_of (deps.begin (), deps.end (), [](const NodeRef& dep) -> bool {
2129 return dep->hasNumericalProperty (NumericalProperty::ConstantZero);
2130 }))
2131 {
2132 return ConstantZero<R>::create (c, dim);
2133 }
2134 // Select node implementation
2135 bool identityDep0 = deps[0]->hasNumericalProperty (NumericalProperty::ConstantIdentity);
2136 bool identityDep1 = deps[1]->hasNumericalProperty (NumericalProperty::ConstantIdentity);
2137 if (identityDep0 && identityDep1)
2138 {
2139 // No specific class for Identity
2140 using namespace numeric;
2141 return NumericConstant<R>::create (c, identity (dim));
2142 }
2143 else if (identityDep0 && !identityDep1)
2144 {
2145 return Convert<R, T1>::create (c, {deps[1]}, dim);
2146 }
2147 else if (!identityDep0 && identityDep1)
2148 {
2149 return Convert<R, T0>::create (c, {deps[0]}, dim);
2150 }
2151 else
2152 {
2153 return cachedAs<Value<R>>(c, std::make_shared<Self>(std::move (deps), dim));
2154 }
2155 }
2156
2158 : Value<R>(std::move (deps)), targetDimension_ (dim)
2159 {}
2160
2161 std::string debugInfo () const override
2162 {
2163 using namespace numeric;
2164 return debug (this->accessValueConst ()) + " targetDim=" + to_string (targetDimension_);
2165 }
2166
2167 std::string shape() const override
2168 {
2169 return "doubleoctagon";
2170 }
2171
2172 std::string color() const override
2173 {
2174 return "#9e9e9e";
2175 }
2176
2177
2178 std::string description() const override
2179 {
2180 return "Matrix Product";
2181 }
2182
2183 // MatrixProduct additional arguments = ().
2184 bool compareAdditionalArguments (const Node_DF& other) const final
2185 {
2186 return dynamic_cast<const Self*>(&other) != nullptr;
2187 }
2188
2189 NodeRef derive (Context& c, const Node_DF& node) final
2190 {
2191 if (&node == this)
2192 {
2194 }
2195 const auto& x0 = this->dependency (0);
2196 const auto& x1 = this->dependency (1);
2197 auto dx0_prod = Self::create (c, {x0->derive (c, node), x1}, targetDimension_);
2198 auto dx1_prod = Self::create (c, {x0, x1->derive (c, node)}, targetDimension_);
2199 return CWiseAdd<R, std::tuple<R, R>>::create (c, {dx0_prod, dx1_prod}, targetDimension_);
2200 }
2201
2203 {
2204 return Self::create (c, std::move (deps), targetDimension_);
2205 }
2206
2207private:
2208 void compute () final
2209 {
2210 auto& result = this->accessValueMutable ();
2211 const auto& x0 = accessValueConstCast<DepT0>(*this->dependency (0));
2212 const auto& x1 = accessValueConstCast<DepT1>(*this->dependency (1));
2213 result.noalias () =
2215#ifdef DEBUG
2216 if ((x1.cols() + x1.rows() < 100 ) && (x0.cols() + x0.rows() < 100))
2217 {
2218 std::cerr << "=== MatrixProd === " << this << std::endl;
2219 std::cerr << "x0= " << x0.rows() << " x " << x0.cols() << std::endl;
2220 std::cerr << "x0= " << NumericalDependencyTransform<T0>::transform (x0) << std::endl;
2221 std::cerr << "x1= " << x1.rows() << " x " << x1.cols() << std::endl;
2222 std::cerr << "x1= " << NumericalDependencyTransform<T1>::transform (x1) << std::endl;
2223 std::cerr << "result= " << result << std::endl;
2224 std::cerr << "=== end MatrixProd === " << this << std::endl << std::endl;
2225 }
2226#endif
2227 }
2228
2230};
2231
2232/*************************************************************************
2233 * @brief r = n * delta + x.
2234 * - r: T.
2235 * - delta: double.
2236 * - x: T.
2237 * - n: constant int.
2238 * - Order of dependencies: (delta, x).
2239 *
2240 * Adds n * delta to all values (component wise) of x.
2241 * Used to generate x +/- delta values for numerical derivation.
2242 * Node construction should be done with the create static method.
2243 */
2244template<typename T> class ShiftDelta : public Value<T>
2245{
2246public:
2248
2250 static ValueRef<T> create (Context& c, NodeRefVec&& deps, int n, const Dimension<T>& dim)
2251 {
2252 // Check dependencies
2253 checkDependenciesNotNull (typeid (Self), deps);
2254 checkDependencyVectorSize (typeid (Self), deps, 2);
2255 checkNthDependencyIsValue<double>(typeid (Self), deps, 0);
2256 checkNthDependencyIsValue<T>(typeid (Self), deps, 1);
2257 // Detect if we have a chain of ShiftDelta with the same delta.
2258 auto& delta = deps[0];
2259 auto& x = deps[1];
2260 auto* xAsShiftDelta = dynamic_cast<const ShiftDelta<T>*>(x.get ());
2261 if (xAsShiftDelta != nullptr && xAsShiftDelta->dependency (0) == delta)
2262 {
2263 // Merge with ShiftDelta dependency by summing the n.
2264 return Self::create (c, NodeRefVec{x->dependencies ()}, n + xAsShiftDelta->getN (), dim);
2265 }
2266 // Not a merge, select node implementation.
2267 if (n == 0 || delta->hasNumericalProperty (NumericalProperty::ConstantZero))
2268 {
2269 return convertRef<Value<T>>(x);
2270 }
2271 else
2272 {
2273 return cachedAs<Value<T>>(c, std::make_shared<Self>(std::move (deps), n, dim));
2274 }
2275 }
2276
2277 ShiftDelta (NodeRefVec&& deps, int n, const Dimension<T>& dim)
2278 : Value<T>(std::move (deps)), targetDimension_ (dim), n_ (n) {}
2279
2280 std::string debugInfo () const override
2281 {
2282 using namespace numeric;
2283 return debug (this->accessValueConst ()) + " targetDim=" + to_string (targetDimension_) +
2284 " n=" + std::to_string (n_);
2285 }
2286
2287 // ShiftDelta additional arguments = (n_).
2288 bool compareAdditionalArguments (const Node_DF& other) const final
2289 {
2290 const auto* derived = dynamic_cast<const Self*>(&other);
2291 return derived != nullptr && n_ == derived->n_;
2292 }
2293 std::size_t hashAdditionalArguments () const final
2294 {
2295 std::size_t seed = 0;
2296 combineHash (seed, n_);
2297 return seed;
2298 }
2299
2300 NodeRef derive (Context& c, const Node_DF& node) final
2301 {
2302 if (&node == this)
2303 {
2305 }
2306 auto& delta = this->dependency (0);
2307 auto& x = this->dependency (1);
2308 return Self::create (c, {delta->derive (c, node), x->derive (c, node)}, n_, targetDimension_);
2309 }
2310
2312 {
2313 return Self::create (c, std::move (deps), n_, targetDimension_);
2314 }
2315
2316 int getN () const { return n_; }
2317
2318private:
2319 void compute () final
2320 {
2321 using namespace numeric;
2322 auto& result = this->accessValueMutable ();
2323 const auto& delta = accessValueConstCast<double>(*this->dependency (0));
2324 const auto& x = accessValueConstCast<T>(*this->dependency (1));
2325 result = n_ * delta + cwise (x);
2326 }
2327
2329 int n_;
2330};
2331
2332
2333/*************************************************************************
2334 * @brief r = (1/delta)^n * sum_i coeffs_i * x_i.
2335 * - r: T.
2336 * - delta: double.
2337 * - x_i: T.
2338 * - n: constant int.
2339 * - coeffs_i: constant double.
2340 * - Order of dependencies: (delta, x_i).
2341 *
2342 * Weighted sum of dependencies, multiplied by a double.
2343 * Used to combine f(x+n*delta) values in numerical derivation.
2344 * Node construction should be done with the create static method.
2345 *
2346 * Lambda represents the 1/delta^n for a nth-order numerical derivation.
2347 * Note that in the whole class, coeffs[i] is the coefficient of deps[i + 1] !
2348 */
2349
2350template<typename T> class CombineDeltaShifted : public Value<T>
2351{
2352public:
2354
2356 static ValueRef<T> create (Context& c, NodeRefVec&& deps, int n, std::vector<double>&& coeffs,
2357 const Dimension<T>& dim)
2358 {
2359 // Check dependencies
2360 checkDependenciesNotNull (typeid (Self), deps);
2361 checkDependencyVectorSize (typeid (Self), deps, 1 + coeffs.size ());
2362 checkNthDependencyIsValue<double>(typeid (Self), deps, 0);
2363 checkDependencyRangeIsValue<T>(typeid (Self), deps, 1, deps.size ());
2364 //
2365 auto cleanAndCreateNode = [&c, &dim](NodeRefVec&& deps2, int n2,
2366 std::vector<double>&& coeffs2) -> ValueRef<T> {
2367 // Clean deps : remove constant 0, of deps with a 0 factor.
2368 for (std::size_t i = 0; i < coeffs2.size ();)
2369 {
2370 if (coeffs2[i] == 0. || deps2[i + 1]->hasNumericalProperty (NumericalProperty::ConstantZero))
2371 {
2372 coeffs2.erase (coeffs2.begin () + std::ptrdiff_t (i));
2373 deps2.erase (deps2.begin () + std::ptrdiff_t (i + 1));
2374 }
2375 else
2376 {
2377 ++i;
2378 }
2379 }
2380 // Final node selection
2381 if (coeffs2.empty ())
2382 {
2383 return ConstantZero<T>::create (c, dim);
2384 }
2385 else
2386 {
2387 return cachedAs<Value<T>>(
2388 c, std::make_shared<Self>(std::move (deps2), n2, std::move (coeffs2), dim));
2389 }
2390 };
2391 // Detect if we can merge this node with its dependencies
2392 const auto& delta = deps[0];
2393 auto isSelfWithSameDelta = [&delta](const NodeRef& dep) -> bool {
2394 return dynamic_cast<const Self*>(dep.get ()) != nullptr && dep->dependency (0) == delta;
2395 };
2396 if (!coeffs.empty () && std::all_of (deps.begin () + 1, deps.end (), isSelfWithSameDelta))
2397 {
2398 const auto depN = static_cast<const Self&>(*deps[1]).getN ();
2399 auto useSameNasDep1 = [depN](const NodeRef& dep) {
2400 return static_cast<const Self&>(*dep).getN () == depN;
2401 };
2402 if (std::all_of (deps.begin () + 2, deps.end (), useSameNasDep1))
2403 {
2404 /* Merge with dependencies because they use the same delta and a common N.
2405 *
2406 * V(this) = 1/delta^n * sum_i V(deps[i]) * coeffs[i].
2407 * V(deps[i]) = 1/delta^depN * sum_j V(depDeps[i][j]) * depCoeffs[i][j].
2408 * V(this) = 1/delta^(n + depN) * sum_ij V(depDeps[i][j]) * coeffs[i] * depCoeffs[i][j].
2409 * And simplify the sum by summing coeffs for each unique depDeps[i][j].
2410 */
2411 NodeRefVec mergedDeps{delta};
2412 std::vector<double> mergedCoeffs;
2413 for (std::size_t i = 0; i < coeffs.size (); ++i)
2414 {
2415 const auto& dep = static_cast<const Self&>(*deps[i + 1]);
2416 const auto& depCoeffs = dep.getCoeffs ();
2417 for (std::size_t j = 0; j < depCoeffs.size (); ++j)
2418 {
2419 const auto& subDep = dep.dependency (j + 1);
2420 auto it = std::find (mergedDeps.begin () + 1, mergedDeps.end (), subDep);
2421 if (it != mergedDeps.end ())
2422 {
2423 // Found
2424 const auto subDepIndexInMerged =
2425 static_cast<std::size_t>(std::distance (mergedDeps.begin () + 1, it));
2426 mergedCoeffs[subDepIndexInMerged] += coeffs[i] * depCoeffs[j];
2427 }
2428 else
2429 {
2430 // Not found
2431 mergedDeps.emplace_back (subDep);
2432 mergedCoeffs.emplace_back (coeffs[i] * depCoeffs[j]);
2433 }
2434 }
2435 }
2436 return cleanAndCreateNode (std::move (mergedDeps), n + depN, std::move (mergedCoeffs));
2437 }
2438 }
2439 // If not able to merge
2440 return cleanAndCreateNode (std::move (deps), n, std::move (coeffs));
2441 }
2442
2443 CombineDeltaShifted (NodeRefVec&& deps, int n, std::vector<double>&& coeffs, const Dimension<T>& dim)
2444 : Value<T>(std::move (deps)), targetDimension_ (dim), coeffs_ (std::move (coeffs)), n_ (n)
2445 {
2446 assert (this->coeffs_.size () + 1 == this->nbDependencies ());
2447 }
2448
2449 std::string description() const override
2450 {
2451 std::string s = " CombineDeltaShifted n=" + std::to_string (n_) + " coeffs={";
2452 if (!coeffs_.empty ())
2453 {
2454 s += std::to_string (coeffs_[0]);
2455 for (std::size_t i = 1; i < coeffs_.size (); ++i)
2456 {
2457 s += ';' + std::to_string (coeffs_[i]);
2458 }
2459 }
2460 s += '}';
2461 return s;
2462 }
2463
2464 std::string debugInfo () const override
2465 {
2466 using namespace numeric;
2467 std::string s = debug (this->accessValueConst ()) + " targetDim=" + to_string (targetDimension_);
2468 return s;
2469 }
2470
2471 // CombineDeltaShifted additional arguments = (n_, coeffs_).
2472 bool compareAdditionalArguments (const Node_DF& other) const final
2473 {
2474 const auto* derived = dynamic_cast<const Self*>(&other);
2475 return derived != nullptr && n_ == derived->n_ && coeffs_ == derived->coeffs_;
2476 }
2477 std::size_t hashAdditionalArguments () const final
2478 {
2479 std::size_t seed = 0;
2480 combineHash (seed, n_);
2481 for (const auto d : coeffs_)
2482 {
2483 combineHash (seed, d);
2484 }
2485 return seed;
2486 }
2487
2488 NodeRef derive (Context& c, const Node_DF& node) final
2489 {
2490 if (&node == this)
2491 {
2493 }
2494 // For simplicity, we assume delta is a constant with respect to derivation node.
2495 auto& delta = this->dependency (0);
2496 if (isTransitivelyDependentOn (node, *delta))
2497 {
2498 // Fail if delta is not constant for node.
2499 failureDeltaNotDerivable (typeid (Self));
2500 }
2501 // Derivation is a simple weighted sums of derivatives.
2502 const auto nbDeps = this->nbDependencies ();
2503 NodeRefVec derivedDeps (nbDeps);
2504 derivedDeps[0] = delta;
2505 for (std::size_t i = 1; i < nbDeps; ++i)
2506 {
2507 derivedDeps[i] = this->dependency (i)->derive (c, node);
2508 }
2509 return Self::create (c, std::move (derivedDeps), n_, std::vector<double>{coeffs_}, targetDimension_);
2510 }
2511
2513 {
2514 return Self::create (c, std::move (deps), n_, std::vector<double>{coeffs_}, targetDimension_);
2515 }
2516
2517 const std::vector<double>& getCoeffs () const { return coeffs_; }
2518
2519 int getN () const { return n_; }
2520
2521private:
2522 void compute() override { compute<T>();}
2523
2524 template<class U>
2525 typename std::enable_if<!std::is_same<U, TransitionFunction>::value, void>::type
2527 {
2528 using namespace numeric;
2529 auto& result = this->accessValueMutable ();
2530 const auto& delta = accessValueConstCast<double>(*this->dependency (0));
2531 const double lambda = pow (delta, -n_);
2532 result = zero (targetDimension_);
2533 for (std::size_t i = 0; i < coeffs_.size (); ++i)
2534 {
2535 const auto& x = accessValueConstCast<T>(*this->dependency (1 + i));
2536 cwise (result) += (lambda * coeffs_[i]) * cwise (x);
2537 }
2538 }
2539
2540 template<class U>
2541 typename std::enable_if<std::is_same<U, TransitionFunction>::value, void>::type
2543 {
2544 using namespace numeric;
2545 auto& result = this->accessValueMutable ();
2546 const auto& delta = accessValueConstCast<double>(*this->dependency (0));
2547 const double lambda = pow (delta, -n_);
2548
2549 std::vector<const T*> vT;
2550 const auto& dep = this->dependencies();
2551
2552 for (size_t i = 1; i < dep.size(); i++)
2553 {
2554 vT.push_back(&accessValueConstCast<T>(*dep[i]));
2555 }
2556
2557 result = [vT, lambda, this](const VectorLik& x)->VectorLik {
2558 using namespace numeric;
2559 VectorLik r = zero (Dimension<VectorLik>(this->targetDimension_.cols, (Eigen::Index)1));
2560
2561 for (std::size_t i = 0; i < this->coeffs_.size (); ++i)
2562 {
2563 cwise(r) += (lambda * this->coeffs_[i]) * cwise((*vT[i])(x));
2564 }
2565 return r;
2566 };
2567 }
2568
2570 std::vector<double> coeffs_;
2571 int n_;
2572};
2573
2574// Precompiled instantiations of numeric nodes
2575
2577
2578extern template class CWiseAdd<double, std::tuple<double, double>>;
2581extern template class CWiseAdd<RowLik, std::tuple<RowLik, RowLik>>;
2584
2585extern template class CWiseAdd<RowLik, MatrixLik>;
2586extern template class CWiseAdd<VectorLik, MatrixLik>;
2587extern template class CWiseAdd<DataLik, VectorLik>;
2588extern template class CWiseAdd<DataLik, RowLik>;
2589
2590extern template class CWiseAdd<double, ReductionOf<double>>;
2592extern template class CWiseAdd<VectorLik, ReductionOf<VectorLik>>;
2593extern template class CWiseAdd<RowLik, ReductionOf<RowLik>>;
2594extern template class CWiseAdd<MatrixLik, ReductionOf<MatrixLik>>;
2596
2602
2603extern template class CWiseMean<VectorLik, ReductionOf<VectorLik>, Eigen::VectorXd>;
2604extern template class CWiseMean<RowLik, ReductionOf<RowLik>, Eigen::VectorXd>;
2605extern template class CWiseMean<MatrixLik, ReductionOf<MatrixLik>, Eigen::VectorXd>;
2606extern template class CWiseMean<VectorLik, ReductionOf<VectorLik>, Eigen::RowVectorXd>;
2607extern template class CWiseMean<RowLik, ReductionOf<RowLik>, Eigen::RowVectorXd>;
2608extern template class CWiseMean<MatrixLik, ReductionOf<MatrixLik>, Eigen::RowVectorXd>;
2609
2610extern template class CWiseSub<double, std::tuple<double, double>>;
2613extern template class CWiseSub<RowLik, std::tuple<RowLik, RowLik>>;
2615
2618
2619extern template class CWiseMul<double, std::tuple<double, double>>;
2620extern template class CWiseMul<double, std::tuple<double, uint>>;
2624extern template class CWiseMul<RowLik, std::tuple<RowLik, RowLik>>;
2626
2634
2635extern template class CWiseMul<double, ReductionOf<double>>;
2637extern template class CWiseMul<VectorLik, ReductionOf<VectorLik>>;
2638extern template class CWiseMul<RowLik, ReductionOf<RowLik>>;
2639extern template class CWiseMul<MatrixLik, ReductionOf<MatrixLik>>;
2640
2641extern template class CWiseNegate<double>;
2642extern template class CWiseNegate<ExtendedFloat>;
2643extern template class CWiseNegate<VectorLik>;
2644extern template class CWiseNegate<RowLik>;
2645extern template class CWiseNegate<MatrixLik>;
2646
2647extern template class CWiseInverse<double>;
2648extern template class CWiseInverse<ExtendedFloat>;
2649extern template class CWiseInverse<VectorLik>;
2650extern template class CWiseInverse<RowLik>;
2651extern template class CWiseInverse<MatrixLik>;
2652
2653extern template class CWiseLog<double>;
2654extern template class CWiseLog<ExtendedFloat>;
2655extern template class CWiseLog<VectorLik>;
2656extern template class CWiseLog<RowLik>;
2657extern template class CWiseLog<MatrixLik>;
2658
2659extern template class CWiseExp<double>;
2660extern template class CWiseExp<ExtendedFloat>;
2661extern template class CWiseExp<VectorLik>;
2662extern template class CWiseExp<RowLik>;
2663extern template class CWiseExp<MatrixLik>;
2664
2665extern template class CWiseConstantPow<double>;
2666extern template class CWiseConstantPow<ExtendedFloat>;
2667extern template class CWiseConstantPow<VectorLik>;
2668extern template class CWiseConstantPow<RowLik>;
2669extern template class CWiseConstantPow<MatrixLik>;
2670
2672extern template class ScalarProduct<DataLik, RowLik, RowLik>;
2673
2676
2677extern template class SumOfLogarithms<VectorLik>;
2678extern template class SumOfLogarithms<RowLik>;
2679
2684
2685extern template class ShiftDelta<double>;
2686extern template class ShiftDelta<VectorLik>;
2687extern template class ShiftDelta<RowLik>;
2688extern template class ShiftDelta<MatrixLik>;
2689
2690extern template class CombineDeltaShifted<double>;
2691extern template class CombineDeltaShifted<VectorLik>;
2692extern template class CombineDeltaShifted<RowLik>;
2693extern template class CombineDeltaShifted<MatrixLik>;
2694extern template class CombineDeltaShifted<TransitionFunction>;
2695
2696
2697/*****************************************************************************
2698 * Numerical derivation helpers.
2699 */
2701
2704{
2707};
2708
2723template<typename NodeT, typename DepT, typename B>
2725 NodeRef dep,
2726 const Dimension<DepT>& depDim,
2727 Dimension<NodeT>& nodeDim,
2728 B buildNodeWithDep)
2729{
2730 if (config.delta == nullptr)
2731 {
2733 }
2734 switch (config.type)
2735 {
2737 // Shift {-1, +1}, coeffs {-0.5, +0.5}
2738 auto shift_m1 = ShiftDelta<DepT>::create (c, {config.delta, dep}, -1, depDim);
2739 auto shift_p1 = ShiftDelta<DepT>::create (c, {config.delta, dep}, 1, depDim);
2740 NodeRefVec combineDeps (3);
2741 combineDeps[0] = config.delta;
2742 combineDeps[1] = buildNodeWithDep (std::move (shift_m1));
2743 combineDeps[2] = buildNodeWithDep (std::move (shift_p1));
2744 return CombineDeltaShifted<NodeT>::create (c, std::move (combineDeps), 1, {-0.5, 0.5}, nodeDim);
2745 } break;
2747 // Shift {-2, -1, +1, +2}, coeffs {1/12, -2/3, 2/3, -1/12}
2748 auto shift_m2 = ShiftDelta<DepT>::create (c, {config.delta, dep}, -2, depDim);
2749 auto shift_m1 = ShiftDelta<DepT>::create (c, {config.delta, dep}, -1, depDim);
2750 auto shift_p1 = ShiftDelta<DepT>::create (c, {config.delta, dep}, 1, depDim);
2751 auto shift_p2 = ShiftDelta<DepT>::create (c, {config.delta, dep}, 2, depDim);
2752 NodeRefVec combineDeps (5);
2753 combineDeps[0] = config.delta;
2754 combineDeps[1] = buildNodeWithDep (std::move (shift_m2));
2755 combineDeps[2] = buildNodeWithDep (std::move (shift_m1));
2756 combineDeps[3] = buildNodeWithDep (std::move (shift_p1));
2757 combineDeps[4] = buildNodeWithDep (std::move (shift_p2));
2758 return CombineDeltaShifted<NodeT>::create (c, std::move (combineDeps), 1,
2759 {1. / 12., -2. / 3., 2. / 3., -1. / 12.}, nodeDim);
2760 } break;
2761 default:
2763 }
2764}
2765
2766template<typename NodeT, typename B>
2769 NodeRef dep,
2770 const Dimension<NodeT>& nodeDim,
2771 B buildNodeWithDep)
2772{
2773 if (config.delta == nullptr)
2774 {
2776 }
2777 ConfiguredParameter* param = dynamic_cast<ConfiguredParameter*>(dep.get());
2778 if (param == nullptr)
2779 throw Exception("generateNumericalDerivative : dependency should be ConfiguredParameter");
2780
2781 switch (config.type)
2782 {
2784 // Shift {-1, +1}, coeffs {-0.5, +0.5}
2785 auto shift_m1 = ShiftParameter::create (c, {dep->dependency(0), config.delta}, *param, -1);
2786 auto shift_p1 = ShiftParameter::create (c, {dep->dependency(0), config.delta}, *param, 1);
2787 NodeRefVec combineDeps (3);
2788 combineDeps[0] = config.delta;
2789 combineDeps[1] = buildNodeWithDep (std::move (shift_m1));
2790 combineDeps[2] = buildNodeWithDep (std::move (shift_p1));
2791 return CombineDeltaShifted<NodeT>::create (c, std::move (combineDeps), 1, {-0.5, 0.5}, nodeDim);
2792 } break;
2794 // Shift {-2, -1, +1, +2}, coeffs {1/12, -2/3, 2/3, -1/12}
2795 auto shift_m2 = ShiftParameter::create (c, {dep->dependency(0), config.delta}, *param, -2);
2796 auto shift_m1 = ShiftParameter::create (c, {dep->dependency(0), config.delta}, *param, -1);
2797 auto shift_p1 = ShiftParameter::create (c, {dep->dependency(0), config.delta}, *param, 1);
2798 auto shift_p2 = ShiftParameter::create (c, {dep->dependency(0), config.delta}, *param, 2);
2799 NodeRefVec combineDeps (5);
2800 combineDeps[0] = config.delta;
2801 combineDeps[1] = buildNodeWithDep (std::move (shift_m2));
2802 combineDeps[2] = buildNodeWithDep (std::move (shift_m1));
2803 combineDeps[3] = buildNodeWithDep (std::move (shift_p1));
2804 combineDeps[4] = buildNodeWithDep (std::move (shift_p2));
2805 return CombineDeltaShifted<NodeT>::create (c, std::move (combineDeps), 1,
2806 {1. / 12., -2. / 3., 2. / 3., -1. / 12.}, nodeDim);
2807 } break;
2808 default:
2810 }
2811}
2812} // namespace bpp
2813#endif // BPP_PHYL_LIKELIHOOD_DATAFLOW_DATAFLOWCWISECOMPUTING_H
CWiseAdd(NodeRefVec &&deps, const Dimension< R > &dim)
std::enable_if< std::is_same< U, TransitionFunction >::value, void >::type compute()
Computation implementation.
bool compareAdditionalArguments(const Node_DF &other) const final
Compare node-specific configuration to another.
std::enable_if<!std::is_same< U, TransitionFunction >::value, void >::type compute()
Computation implementation.
std::string debugInfo() const override
Node debug info (default = ""): user defined detailed info for DF graph debug.
NodeRef recreate(Context &c, NodeRefVec &&deps) final
Recreate the node with different dependencies.
NodeRef derive(Context &c, const Node_DF &node) final
Returns a node computing d(this_node_expression)/d(node_expression).
static ValueRef< R > create(Context &c, NodeRefVec &&deps, const Dimension< R > &dim)
Build a new CWiseAdd node with the given output dimensions.
void compute() override
Computation implementation.
void compute() override
Computation implementation.
bool compareAdditionalArguments(const Node_DF &other) const final
Compare node-specific configuration to another.
NodeRef derive(Context &c, const Node_DF &node) final
Returns a node computing d(this_node_expression)/d(node_expression).
std::string debugInfo() const override
Node debug info (default = ""): user defined detailed info for DF graph debug.
std::enable_if<!std::is_same< U, TransitionFunction >::value, void >::type compute()
Computation implementation.
std::string description() const override
Node pretty name (default = type name).
static ValueRef< R > create(Context &c, NodeRefVec &&deps, const Dimension< R > &dim)
Build a new CWiseAdd node with the given output dimensions.
NodeRef recreate(Context &c, NodeRefVec &&deps) final
Recreate the node with different dependencies.
CWiseAdd(NodeRefVec &&deps, const Dimension< R > &dim)
std::enable_if< std::is_same< U, TransitionFunction >::value, void >::type compute()
Computation implementation.
std::enable_if<(std::is_base_of< V, MatrixLik >::value)&&(std::is_same< U, VectorLik >::value||std::is_same< U, Eigen::VectorXd >::value), void >::type compute()
Computation implementation.
NodeRef derive(Context &c, const Node_DF &node) final
Returns a node computing d(this_node_expression)/d(node_expression).
std::enable_if<(std::is_base_of< V, VectorLik >::value)||(std::is_same< V, RowLik >::value||std::is_same< V, Eigen::RowVectorXd >::value), void >::type compute()
Computation implementation.
void compute() override
Computation implementation.
std::enable_if<(std::is_base_of< V, MatrixLik >::value)&&(std::is_same< U, RowLik >::value||std::is_same< U, Eigen::RowVectorXd >::value), void >::type compute()
Computation implementation.
std::string debugInfo() const override
Node debug info (default = ""): user defined detailed info for DF graph debug.
NodeRef recreate(Context &c, NodeRefVec &&deps) final
Recreate the node with different dependencies.
Dimension< R > targetDimension_
static ValueRef< R > create(Context &c, NodeRefVec &&deps, const Dimension< R > &dim)
Build a new CWiseAdd node with the given output dimensions.
CWiseAdd(NodeRefVec &&deps, const Dimension< R > &dim)
bool compareAdditionalArguments(const Node_DF &other) const final
Compare node-specific configuration to another.
void compute() override
Computation implementation.
std::enable_if< std::is_base_of< U, MatrixLik >::value &&std::is_same< F, TransitionFunction >::value, void >::type compute()
Computation implementation.
bool compareAdditionalArguments(const Node_DF &other) const final
Compare node-specific configuration to another.
static ValueRef< R > create(Context &c, NodeRefVec &&deps, const Dimension< R > &dim)
Build a new CWiseApply node.
std::string debugInfo() const override
Node debug info (default = ""): user defined detailed info for DF graph debug.
NodeRef recreate(Context &c, NodeRefVec &&deps) final
Recreate the node with different dependencies.
std::string shape() const override
std::string color() const override
CWiseApply(NodeRefVec &&deps, const Dimension< R > &dim)
Dimension< R > targetDimension_
std::vector< VectorLik > bppLik_
For computation purpose.
NodeRef derive(Context &c, const Node_DF &node) final
Returns a node computing d(this_node_expression)/d(node_expression).
std::string description() const override
Node pretty name (default = type name).
NodeRef derive(Context &c, const Node_DF &node) final
Returns a node computing d(this_node_expression)/d(node_expression).
bool compareAdditionalArguments(const Node_DF &other) const final
Compare node-specific configuration to another.
static ValueRef< T > create(Context &c, NodeRefVec &&deps, double exponent, double factor, const Dimension< T > &dim)
Build a new CWiseConstantPow node with the given output dimensions and factors.
void compute() final
Computation implementation.
NodeRef recreate(Context &c, NodeRefVec &&deps) final
Recreate the node with different dependencies.
std::size_t hashAdditionalArguments() const final
Return the hash of node-specific configuration.
std::string debugInfo() const override
Node debug info (default = ""): user defined detailed info for DF graph debug.
CWiseConstantPow(NodeRefVec &&deps, double exponent, double factor, const Dimension< T > &dim)
std::enable_if< std::is_same< U, TransitionFunction >::value &&std::is_same< V, TransitionFunction >::value, void >::type compute()
Computation implementation.
static ValueRef< R > create(Context &c, NodeRefVec &&deps, const Dimension< R > &dim)
Build a new CWiseDiv node with the given output dimensions.
NodeRef derive(Context &c, const Node_DF &node) final
Returns a node computing d(this_node_expression)/d(node_expression).
std::enable_if< std::is_same< U, TransitionFunction >::value &&!std::is_same< V, TransitionFunction >::value, void >::type compute()
Computation implementation.
std::enable_if<!std::is_same< U, TransitionFunction >::value &&std::is_same< V, TransitionFunction >::value, void >::type compute()
Computation implementation.
std::enable_if<!std::is_same< U, TransitionFunction >::value &&!std::is_same< V, TransitionFunction >::value, void >::type compute()
Computation implementation.
std::string debugInfo() const override
Node debug info (default = ""): user defined detailed info for DF graph debug.
std::string description() const override
Node pretty name (default = type name).
NodeRef recreate(Context &c, NodeRefVec &&deps) final
Recreate the node with different dependencies.
CWiseDiv(NodeRefVec &&deps, const Dimension< R > &dim)
void compute() override
Computation implementation.
bool compareAdditionalArguments(const Node_DF &other) const final
Compare node-specific configuration to another.
Dimension< T > targetDimension_
std::string debugInfo() const override
Node debug info (default = ""): user defined detailed info for DF graph debug.
void compute() final
Computation implementation.
NodeRef derive(Context &c, const Node_DF &node) final
Returns a node computing d(this_node_expression)/d(node_expression).
static ValueRef< T > create(Context &c, NodeRefVec &&deps, const Dimension< T > &dim)
Build a new CWiseExp node with the given output dimensions.
NodeRef recreate(Context &c, NodeRefVec &&deps) final
Recreate the node with different dependencies.
bool compareAdditionalArguments(const Node_DF &other) const final
Compare node-specific configuration to another.
CWiseExp(NodeRefVec &&deps, const Dimension< T > &dim)
CWiseInverse(NodeRefVec &&deps, const Dimension< T > &dim)
NodeRef recreate(Context &c, NodeRefVec &&deps) final
Recreate the node with different dependencies.
bool compareAdditionalArguments(const Node_DF &other) const final
Compare node-specific configuration to another.
NodeRef derive(Context &c, const Node_DF &node) final
Returns a node computing d(this_node_expression)/d(node_expression).
Dimension< T > targetDimension_
static ValueRef< T > create(Context &c, NodeRefVec &&deps, const Dimension< T > &dim)
Build a new CWiseInverse node with the given output dimensions.
void compute() final
Computation implementation.
std::string debugInfo() const override
Node debug info (default = ""): user defined detailed info for DF graph debug.
NodeRef derive(Context &c, const Node_DF &node) final
Returns a node computing d(this_node_expression)/d(node_expression).
CWiseLog(NodeRefVec &&deps, const Dimension< T > &dim)
std::string debugInfo() const override
Node debug info (default = ""): user defined detailed info for DF graph debug.
bool compareAdditionalArguments(const Node_DF &other) const final
Compare node-specific configuration to another.
NodeRef recreate(Context &c, NodeRefVec &&deps) final
Recreate the node with different dependencies.
static ValueRef< T > create(Context &c, NodeRefVec &&deps, const Dimension< T > &dim)
Build a new CWiseLog node with the given output dimensions.
void compute() final
Computation implementation.
Dimension< T > targetDimension_
static ValueRef< R > create(Context &c, NodeRefVec &&deps, const Dimension< R > &dim)
Build a new CWiseMean node with the given output dimensions.
std::string debugInfo() const override
Node debug info (default = ""): user defined detailed info for DF graph debug.
void compute() final
Computation implementation.
CWiseMean(NodeRefVec &&deps, const Dimension< R > &dim)
NodeRef derive(Context &c, const Node_DF &node) final
Returns a node computing d(this_node_expression)/d(node_expression).
bool compareAdditionalArguments(const Node_DF &other) const final
Compare node-specific configuration to another.
NodeRef recreate(Context &c, NodeRefVec &&deps) final
Recreate the node with different dependencies.
std::string description() const override
Node pretty name (default = type name).
bool compareAdditionalArguments(const Node_DF &other) const final
Compare node-specific configuration to another.
NodeRef recreate(Context &c, NodeRefVec &&deps) final
Recreate the node with different dependencies.
std::string description() const override
Node pretty name (default = type name).
CWiseMean(NodeRefVec &&deps, const Dimension< R > &dim)
std::string debugInfo() const override
Node debug info (default = ""): user defined detailed info for DF graph debug.
static ValueRef< R > create(Context &c, NodeRefVec &&deps, const Dimension< R > &dim)
Build a new CWiseMean node with the given output dimensions.
NodeRef derive(Context &c, const Node_DF &node) final
Returns a node computing d(this_node_expression)/d(node_expression).
static ValueRef< R > create(Context &c, NodeRefVec &&deps, const Dimension< R > &dim)
Build a new CWiseMul node with the given output dimensions.
std::string debugInfo() const override
Node debug info (default = ""): user defined detailed info for DF graph debug.
bool compareAdditionalArguments(const Node_DF &other) const final
Compare node-specific configuration to another.
NodeRef derive(Context &c, const Node_DF &node) final
Returns a node computing d(this_node_expression)/d(node_expression).
CWiseMul(NodeRefVec &&deps, const Dimension< R > &dim)
NodeRef recreate(Context &c, NodeRefVec &&deps) final
Recreate the node with different dependencies.
void compute() final
Computation implementation.
std::enable_if< std::is_same< U, TransitionFunction >::value &&std::is_same< V, TransitionFunction >::value, void >::type compute()
Computation implementation.
bool compareAdditionalArguments(const Node_DF &other) const final
Compare node-specific configuration to another.
NodeRef recreate(Context &c, NodeRefVec &&deps) final
Recreate the node with different dependencies.
std::enable_if<!std::is_same< U, TransitionFunction >::value &&!std::is_same< V, TransitionFunction >::value, void >::type compute()
Computation implementation.
CWiseMul(NodeRefVec &&deps, const Dimension< R > &dim)
std::string debugInfo() const override
Node debug info (default = ""): user defined detailed info for DF graph debug.
void compute() override
Computation implementation.
static ValueRef< R > create(Context &c, NodeRefVec &&deps, const Dimension< R > &dim)
Build a new CWiseMul node with the given output dimensions.
std::string description() const override
Node pretty name (default = type name).
std::enable_if<!std::is_same< U, TransitionFunction >::value &&std::is_same< V, TransitionFunction >::value, void >::type compute()
Computation implementation.
std::enable_if< std::is_same< U, TransitionFunction >::value &&!std::is_same< V, TransitionFunction >::value, void >::type compute()
Computation implementation.
NodeRef derive(Context &c, const Node_DF &node) final
Returns a node computing d(this_node_expression)/d(node_expression).
CWiseNegate(NodeRefVec &&deps, const Dimension< T > &dim)
static ValueRef< T > create(Context &c, NodeRefVec &&deps, const Dimension< T > &dim)
Build a new CWiseNegate node with the given output dimensions.
bool compareAdditionalArguments(const Node_DF &other) const final
Compare node-specific configuration to another.
void compute() final
Computation implementation.
NodeRef derive(Context &c, const Node_DF &node) final
Returns a node computing d(this_node_expression)/d(node_expression).
NodeRef recreate(Context &c, NodeRefVec &&deps) final
Recreate the node with different dependencies.
Dimension< T > targetDimension_
std::string debugInfo() const override
Node debug info (default = ""): user defined detailed info for DF graph debug.
std::string debugInfo() const override
Node debug info (default = ""): user defined detailed info for DF graph debug.
NodeRef derive(Context &c, const Node_DF &node) final
Returns a node computing d(this_node_expression)/d(node_expression).
void compute() final
Computation implementation.
bool compareAdditionalArguments(const Node_DF &other) const final
Compare node-specific configuration to another.
CWiseSub(NodeRefVec &&deps, const Dimension< R > &dim)
NodeRef recreate(Context &c, NodeRefVec &&deps) final
Recreate the node with different dependencies.
static ValueRef< R > create(Context &c, NodeRefVec &&deps, const Dimension< R > &dim)
Build a new CWiseSub node with the given output dimensions.
std::string debugInfo() const override
Node debug info (default = ""): user defined detailed info for DF graph debug.
std::size_t hashAdditionalArguments() const final
Return the hash of node-specific configuration.
void compute() override
Computation implementation.
static ValueRef< T > create(Context &c, NodeRefVec &&deps, int n, std::vector< double > &&coeffs, const Dimension< T > &dim)
Build a new CombineDeltaShifted node with the given output dimensions, exponent and weights.
std::string description() const override
Node pretty name (default = type name).
std::enable_if<!std::is_same< U, TransitionFunction >::value, void >::type compute()
Computation implementation.
const std::vector< double > & getCoeffs() const
NodeRef recreate(Context &c, NodeRefVec &&deps) final
Recreate the node with different dependencies.
NodeRef derive(Context &c, const Node_DF &node) final
Returns a node computing d(this_node_expression)/d(node_expression).
std::enable_if< std::is_same< U, TransitionFunction >::value, void >::type compute()
Computation implementation.
CombineDeltaShifted(NodeRefVec &&deps, int n, std::vector< double > &&coeffs, const Dimension< T > &dim)
bool compareAdditionalArguments(const Node_DF &other) const final
Compare node-specific configuration to another.
Data flow node representing a parameter.
Definition: Parameter.h:27
static std::shared_ptr< Self > create(Context &c, const Dimension< T > &dim)
Build a new ConstantOne node of the given dimension.
static std::shared_ptr< Self > create(Context &c, const Dimension< T > &dim)
Build a new ConstantZero node of the given dimension.
Context for dataflow node construction.
Definition: DataFlow.h:527
static ValueRef< R > create(Context &c, NodeRefVec &&deps, const Dimension< R > &dim)
Build a new Convert node with the given output dimensions.
double log() const
const ExtType & exponent_part() const noexcept
Definition: ExtendedFloat.h:89
void normalize() noexcept
static const double ln_radix
Definition: ExtendedFloat.h:51
const FloatType & float_part() const noexcept
Definition: ExtendedFloat.h:88
static ExtendedFloat denorm_mul(const ExtendedFloat &lhs, const ExtendedFloat &rhs)
static ValueRef< R > create(Context &c, NodeRefVec &&deps, const Dimension< T0 > &mDim)
Build a new LogSumExp node with the given input matrix dimensions.
std::string debugInfo() const override
Node debug info (default = ""): user defined detailed info for DF graph debug.
NodeRef derive(Context &c, const Node_DF &node) final
Returns a node computing d(this_node_expression)/d(node_expression).
bool compareAdditionalArguments(const Node_DF &other) const final
Compare node-specific configuration to another.
LogSumExp(NodeRefVec &&deps, const Dimension< T0 > &mDim)
void compute() final
Computation implementation.
Dimension< T0 > mTargetDimension_
NodeRef recreate(Context &c, NodeRefVec &&deps) final
Recreate the node with different dependencies.
std::string description() const override
Node pretty name (default = type name).
bool compareAdditionalArguments(const Node_DF &other) const final
Compare node-specific configuration to another.
NodeRef recreate(Context &c, NodeRefVec &&deps) final
Recreate the node with different dependencies.
std::string shape() const override
MatrixProduct(NodeRefVec &&deps, const Dimension< R > &dim)
typename NumericalDependencyTransform< T0 >::DepType DepT0
void compute() final
Computation implementation.
NodeRef derive(Context &c, const Node_DF &node) final
Returns a node computing d(this_node_expression)/d(node_expression).
typename NumericalDependencyTransform< T1 >::DepType DepT1
static ValueRef< R > create(Context &c, NodeRefVec &&deps, const Dimension< R > &dim)
Build a new MatrixProduct node with the given output dimensions.
std::string debugInfo() const override
Node debug info (default = ""): user defined detailed info for DF graph debug.
std::string color() const override
Base dataflow Node class.
Definition: DataFlow.h:152
const NodeRefVec & dependencies() const noexcept
Definition: DataFlow.h:183
const NodeRef & dependency(std::size_t i) const noexcept
Definition: DataFlow.h:185
virtual bool hasNumericalProperty(NumericalProperty prop) const
Test if the node has the given numerical property.
Definition: DataFlow.cpp:168
std::size_t nbDependencies() const noexcept
Number of dependencies (ie nodes we depend on)
Definition: DataFlow.h:181
static std::shared_ptr< Self > create(Context &c, Args &&... args)
Build a new NumericConstant node with T(args...) value.
NodeRef derive(Context &c, const Node_DF &node) final
Returns a node computing d(this_node_expression)/d(node_expression).
void compute() final
Computation implementation.
ScalarProduct(NodeRefVec &&deps)
NodeRef recreate(Context &c, NodeRefVec &&deps) final
Recreate the node with different dependencies.
bool compareAdditionalArguments(const Node_DF &other) const final
Compare node-specific configuration to another.
std::string debugInfo() const override
Node debug info (default = ""): user defined detailed info for DF graph debug.
static ValueRef< R > create(Context &c, NodeRefVec &&deps)
Build a new ScalarProduct node.
Dimension< T > targetDimension_
void compute() final
Computation implementation.
NodeRef derive(Context &c, const Node_DF &node) final
Returns a node computing d(this_node_expression)/d(node_expression).
bool compareAdditionalArguments(const Node_DF &other) const final
Compare node-specific configuration to another.
std::string debugInfo() const override
Node debug info (default = ""): user defined detailed info for DF graph debug.
ShiftDelta(NodeRefVec &&deps, int n, const Dimension< T > &dim)
std::size_t hashAdditionalArguments() const final
Return the hash of node-specific configuration.
static ValueRef< T > create(Context &c, NodeRefVec &&deps, int n, const Dimension< T > &dim)
Build a new ShiftDelta node with the given output dimensions and shift number.
NodeRef recreate(Context &c, NodeRefVec &&deps) final
Recreate the node with different dependencies.
static std::shared_ptr< ConfiguredParameter > create(Context &c, NodeRefVec &&deps, Parameter &param, const int n)
Build a new ShiftDelta node with the given output dimensions and shift number.
Definition: Parameter.h:136
bool compareAdditionalArguments(const Node_DF &other) const final
Compare node-specific configuration to another.
SumOfLogarithms(NodeRefVec &&deps, const Dimension< F > &mDim)
std::enable_if< std::is_convertible< G *, Eigen::DenseBase< G > * >::value, void >::type compute()
Computation implementation.
NodeRef derive(Context &c, const Node_DF &node) final
Returns a node computing d(this_node_expression)/d(node_expression).
std::enable_if< std::is_convertible< G *, ExtendedFloatEigenBase< G > * >::value, void >::type compute()
Computation implementation.
static ValueRef< DataLik > create(Context &c, NodeRefVec &&deps, const Dimension< F > &mDim)
Build a new SumOfLogarithms node with the given input matrix dimensions.
std::string debugInfo() const override
Node debug info (default = ""): user defined detailed info for DF graph debug.
void compute() override
Computation implementation.
NodeRef recreate(Context &c, NodeRefVec &&deps) final
Recreate the node with different dependencies.
Eigen::Matrix< ExtendedFloat, 1, Eigen::Dynamic > temp_
Abstract Node storing a value of type T.
Definition: DataFlow.h:352
R & accessValueMutable() noexcept
Definition: DataFlow.h:416
const R & accessValueConst() const noexcept
Raw value access (const).
Definition: DataFlow.h:385
T identity(const Dimension< T > &)
bool isinf(const double &d)
T one(const Dimension< T > &)
T & cwise(T &t)
Definition: DataFlowCWise.h:29
std::string debug(const T &t, typename std::enable_if< std::is_arithmetic< T >::value >::type *=0)
T zero(const Dimension< T > &)
Defines the basic types of data flow nodes.
void removeDependenciesIf(NodeRefVec &deps, Predicate p)
double log(const ExtendedFloat &ef)
ExtendedFloat exp(const ExtendedFloat &ef)
NumericConstant< char > NodeX('X')
std::shared_ptr< Value< T > > ValueRef
Shared pointer alias for Value<T>.
Definition: DataFlow.h:84
void checkDependenciesNotNull(const std::type_info &contextNodeType, const NodeRefVec &deps)
Checks that all dependencies are not null, throws if not.
Definition: DataFlow.cpp:103
std::string to_string(const NoDimension &)
ExtendedFloat pow(const ExtendedFloat &ef, double exp)
void failureDeltaNotDerivable(const std::type_info &contextNodeType)
std::vector< NodeRef > NodeRefVec
Alias for a dependency vector (of NodeRef).
Definition: DataFlow.h:81
template void copyBppToEigen(const std::vector< ExtendedFloat > &bppVector, Eigen::RowVectorXd &eigenVector)
ValueRef< NodeT > generateNumericalDerivative(Context &c, const NumericalDerivativeConfiguration &config, NodeRef dep, const Dimension< DepT > &depDim, Dimension< NodeT > &nodeDim, B buildNodeWithDep)
Helper used to generate data flow expressions computing a numerical derivative.
void checkDependencyVectorSize(const std::type_info &contextNodeType, const NodeRefVec &deps, std::size_t expectedSize)
Definition: DataFlow.cpp:83
bool isTransitivelyDependentOn(const Node_DF &searchedDependency, const Node_DF &node)
Check if searchedDependency if a transitive dependency of node.
Definition: DataFlow.cpp:254
void checkDependencyVectorMinSize(const std::type_info &contextNodeType, const NodeRefVec &deps, std::size_t expectedMinSize)
Checks the minimum size of a dependency vector, throws if mismatch.
Definition: DataFlow.cpp:93
void combineHash(std::size_t &seed, const T &t)
Combine hashable value to a hash, from Boost library.
Definition: DataFlow.h:33
void failureNumericalDerivationNotConfigured()
std::shared_ptr< Node_DF > NodeRef
Definition: DataFlow.h:78
static const DepType & transform(const DepType &d)
Transform the DepType value: do nothing in the default case.
T DepType
Real type of the dependency: T in the default case.
Configuration for a numerical derivation: what delta to use, and type of derivation.