21 #include <boost/function.hpp> 23 #pragma GCC diagnostic push 24 #pragma GCC diagnostic ignored "-Wunused-variable" 26 #include <boost/bind.hpp> 28 #pragma GCC diagnostic pop 63 template<
class Y,
class X=
double>
77 template <class X, int N = traits<X>::dimension>
79 boost::function<
double(
const X&)> h,
const X& x,
double delta = 1e-5) {
80 double factor = 1.0 / (2.0 * delta);
82 BOOST_STATIC_ASSERT_MSG(
84 "Template argument X must be a manifold type.");
85 BOOST_STATIC_ASSERT_MSG(N>0,
"Template argument X must be fixed-size type or N must be specified.");
91 Eigen::Matrix<double,N,1> g;
93 for (
int j = 0; j < N; j++) {
99 g(j) = (hxplus - hxmin) * factor;
116 template <class Y, class X, int N = traits<X>::dimension>
119 boost::function<Y(
const X&)> h,
const X& x,
double delta = 1e-5) {
120 typedef typename internal::FixedSizeMatrix<Y,X>::type Matrix;
123 "Template argument Y must be a manifold type.");
127 "Template argument X must be a manifold type.");
128 BOOST_STATIC_ASSERT_MSG(N>0,
"Template argument X must be fixed-size type or N must be specified.");
135 const typename TraitsY::TangentVector zeroY = TraitsY::Local(hx, hx);
136 const size_t m = zeroY.size();
139 Eigen::Matrix<double, N, 1> dx;
143 Matrix H = Matrix::Zero(m, N);
144 const double factor = 1.0 / (2.0 * delta);
145 for (
int j = 0; j < N; j++) {
147 const auto dy1 = TraitsY::Local(hx, h(TraitsX::Retract(x, dx)));
149 const auto dy2 = TraitsY::Local(hx, h(TraitsX::Retract(x, dx)));
151 H.col(j) << (dy1 - dy2) * factor;
157 template<
class Y,
class X>
159 double delta = 1e-5) {
160 return numericalDerivative11<Y, X>(boost::bind(h, _1), x, delta);
171 template<
class Y,
class X1,
class X2>
172 typename internal::FixedSizeMatrix<Y,X1>::type
numericalDerivative21(
const boost::function<Y(
const X1&,
const X2&)>& h,
173 const X1& x1,
const X2& x2,
double delta = 1e-5) {
175 "Template argument Y must be a manifold type.");
177 "Template argument X1 must be a manifold type.");
178 return numericalDerivative11<Y, X1>(boost::bind(h, _1, boost::cref(x2)), x1, delta);
182 template<
class Y,
class X1,
class X2>
184 const X2& x2,
double delta = 1e-5) {
185 return numericalDerivative21<Y, X1, X2>(boost::bind(h, _1, _2), x1, x2, delta);
196 template<
class Y,
class X1,
class X2>
198 const X1& x1,
const X2& x2,
double delta = 1e-5) {
202 "Template argument X2 must be a manifold type.");
203 return numericalDerivative11<Y, X2>(boost::bind(h, boost::cref(x1), _1), x2, delta);
207 template<
class Y,
class X1,
class X2>
209 const X2& x2,
double delta = 1e-5) {
210 return numericalDerivative22<Y, X1, X2>(boost::bind(h, _1, _2), x1, x2, delta);
223 template<
class Y,
class X1,
class X2,
class X3>
225 boost::function<Y(
const X1&,
const X2&,
const X3&)> h,
const X1& x1,
226 const X2& x2,
const X3& x3,
double delta = 1e-5) {
228 "Template argument Y must be a manifold type.");
230 "Template argument X1 must be a manifold type.");
231 return numericalDerivative11<Y, X1>(boost::bind(h, _1, boost::cref(x2), boost::cref(x3)), x1, delta);
234 template<
class Y,
class X1,
class X2,
class X3>
235 typename internal::FixedSizeMatrix<Y,X1>::type
numericalDerivative31(Y (*h)(
const X1&,
const X2&,
const X3&),
236 const X1& x1,
const X2& x2,
const X3& x3,
double delta = 1e-5) {
237 return numericalDerivative31<Y, X1, X2, X3>(boost::bind(h, _1, _2, _3), x1,
251 template<
class Y,
class X1,
class X2,
class X3>
253 boost::function<Y(
const X1&,
const X2&,
const X3&)> h,
const X1& x1,
254 const X2& x2,
const X3& x3,
double delta = 1e-5) {
256 "Template argument Y must be a manifold type.");
258 "Template argument X2 must be a manifold type.");
259 return numericalDerivative11<Y, X2>(boost::bind(h, boost::cref(x1), _1, boost::cref(x3)), x2, delta);
262 template<
class Y,
class X1,
class X2,
class X3>
263 inline typename internal::FixedSizeMatrix<Y,X2>::type
numericalDerivative32(Y (*h)(
const X1&,
const X2&,
const X3&),
264 const X1& x1,
const X2& x2,
const X3& x3,
double delta = 1e-5) {
265 return numericalDerivative32<Y, X1, X2, X3>(boost::bind(h, _1, _2, _3), x1,
279 template<
class Y,
class X1,
class X2,
class X3>
281 boost::function<Y(
const X1&,
const X2&,
const X3&)> h,
const X1& x1,
282 const X2& x2,
const X3& x3,
double delta = 1e-5) {
284 "Template argument Y must be a manifold type.");
286 "Template argument X3 must be a manifold type.");
287 return numericalDerivative11<Y, X3>(boost::bind(h, boost::cref(x1), boost::cref(x2), _1), x3, delta);
290 template<
class Y,
class X1,
class X2,
class X3>
291 inline typename internal::FixedSizeMatrix<Y,X3>::type
numericalDerivative33(Y (*h)(
const X1&,
const X2&,
const X3&),
292 const X1& x1,
const X2& x2,
const X3& x3,
double delta = 1e-5) {
293 return numericalDerivative33<Y, X1, X2, X3>(boost::bind(h, _1, _2, _3), x1,
307 template<
class Y,
class X1,
class X2,
class X3,
class X4>
309 boost::function<Y(
const X1&,
const X2&,
const X3&,
const X4&)> h,
const X1& x1,
310 const X2& x2,
const X3& x3,
const X4& x4,
double delta = 1e-5) {
312 "Template argument Y must be a manifold type.");
314 "Template argument X1 must be a manifold type.");
315 return numericalDerivative11<Y, X1>(boost::bind(h, _1, boost::cref(x2), boost::cref(x3), boost::cref(x4)), x1, delta);
318 template<
class Y,
class X1,
class X2,
class X3,
class X4>
319 inline typename internal::FixedSizeMatrix<Y,X1>::type
numericalDerivative41(Y (*h)(
const X1&,
const X2&,
const X3&,
const X4&),
320 const X1& x1,
const X2& x2,
const X3& x3,
const X4& x4,
double delta = 1e-5) {
321 return numericalDerivative41<Y, X1, X2, X3, X4>(boost::bind(h, _1, _2, _3, _4), x1, x2, x3, x4);
334 template<
class Y,
class X1,
class X2,
class X3,
class X4>
336 boost::function<Y(
const X1&,
const X2&,
const X3&,
const X4&)> h,
const X1& x1,
337 const X2& x2,
const X3& x3,
const X4& x4,
double delta = 1e-5) {
339 "Template argument Y must be a manifold type.");
341 "Template argument X2 must be a manifold type.");
342 return numericalDerivative11<Y, X2>(boost::bind(h, boost::cref(x1), _1, boost::cref(x3), boost::cref(x4)), x2, delta);
345 template<
class Y,
class X1,
class X2,
class X3,
class X4>
346 inline typename internal::FixedSizeMatrix<Y,X2>::type
numericalDerivative42(Y (*h)(
const X1&,
const X2&,
const X3&,
const X4&),
347 const X1& x1,
const X2& x2,
const X3& x3,
const X4& x4,
double delta = 1e-5) {
348 return numericalDerivative42<Y, X1, X2, X3, X4>(boost::bind(h, _1, _2, _3, _4), x1, x2, x3, x4);
361 template<
class Y,
class X1,
class X2,
class X3,
class X4>
363 boost::function<Y(
const X1&,
const X2&,
const X3&,
const X4&)> h,
const X1& x1,
364 const X2& x2,
const X3& x3,
const X4& x4,
double delta = 1e-5) {
366 "Template argument Y must be a manifold type.");
368 "Template argument X3 must be a manifold type.");
369 return numericalDerivative11<Y, X3>(boost::bind(h, boost::cref(x1), boost::cref(x2), _1, boost::cref(x4)), x3, delta);
372 template<
class Y,
class X1,
class X2,
class X3,
class X4>
373 inline typename internal::FixedSizeMatrix<Y,X3>::type
numericalDerivative43(Y (*h)(
const X1&,
const X2&,
const X3&,
const X4&),
374 const X1& x1,
const X2& x2,
const X3& x3,
const X4& x4,
double delta = 1e-5) {
375 return numericalDerivative43<Y, X1, X2, X3, X4>(boost::bind(h, _1, _2, _3, _4), x1, x2, x3, x4);
388 template<
class Y,
class X1,
class X2,
class X3,
class X4>
390 boost::function<Y(
const X1&,
const X2&,
const X3&,
const X4&)> h,
const X1& x1,
391 const X2& x2,
const X3& x3,
const X4& x4,
double delta = 1e-5) {
393 "Template argument Y must be a manifold type.");
395 "Template argument X4 must be a manifold type.");
396 return numericalDerivative11<Y, X4>(boost::bind(h, boost::cref(x1), boost::cref(x2), boost::cref(x3), _1), x4, delta);
399 template<
class Y,
class X1,
class X2,
class X3,
class X4>
400 inline typename internal::FixedSizeMatrix<Y,X4>::type
numericalDerivative44(Y (*h)(
const X1&,
const X2&,
const X3&,
const X4&),
401 const X1& x1,
const X2& x2,
const X3& x3,
const X4& x4,
double delta = 1e-5) {
402 return numericalDerivative44<Y, X1, X2, X3, X4>(boost::bind(h, _1, _2, _3, _4), x1, x2, x3, x4);
416 template<
class Y,
class X1,
class X2,
class X3,
class X4,
class X5>
418 boost::function<Y(
const X1&,
const X2&,
const X3&,
const X4&,
const X5&)> h,
const X1& x1,
419 const X2& x2,
const X3& x3,
const X4& x4,
const X5& x5,
double delta = 1e-5) {
421 "Template argument Y must be a manifold type.");
423 "Template argument X1 must be a manifold type.");
424 return numericalDerivative11<Y, X1>(boost::bind(h, _1, boost::cref(x2), boost::cref(x3), boost::cref(x4), boost::cref(x5)), x1, delta);
427 template<
class Y,
class X1,
class X2,
class X3,
class X4,
class X5>
428 inline typename internal::FixedSizeMatrix<Y,X1>::type
numericalDerivative51(Y (*h)(
const X1&,
const X2&,
const X3&,
const X4&,
const X5&),
429 const X1& x1,
const X2& x2,
const X3& x3,
const X4& x4,
const X5& x5,
double delta = 1e-5) {
430 return numericalDerivative51<Y, X1, X2, X3, X4, X5>(boost::bind(h, _1, _2, _3, _4, _5), x1, x2, x3, x4, x5);
444 template<
class Y,
class X1,
class X2,
class X3,
class X4,
class X5>
446 boost::function<Y(
const X1&,
const X2&,
const X3&,
const X4&,
const X5&)> h,
const X1& x1,
447 const X2& x2,
const X3& x3,
const X4& x4,
const X5& x5,
double delta = 1e-5) {
449 "Template argument Y must be a manifold type.");
451 "Template argument X1 must be a manifold type.");
452 return numericalDerivative11<Y, X2>(boost::bind(h, boost::cref(x1), _1, boost::cref(x3), boost::cref(x4), boost::cref(x5)), x2, delta);
455 template<
class Y,
class X1,
class X2,
class X3,
class X4,
class X5>
456 inline typename internal::FixedSizeMatrix<Y,X2>::type
numericalDerivative51(Y (*h)(
const X1&,
const X2&,
const X3&,
const X4&,
const X5&),
457 const X1& x1,
const X2& x2,
const X3& x3,
const X4& x4,
const X5& x5,
double delta = 1e-5) {
458 return numericalDerivative52<Y, X1, X2, X3, X4, X5>(boost::bind(h, _1, _2, _3, _4, _5), x1, x2, x3, x4, x5);
472 template<
class Y,
class X1,
class X2,
class X3,
class X4,
class X5>
474 boost::function<Y(
const X1&,
const X2&,
const X3&,
const X4&,
const X5&)> h,
const X1& x1,
475 const X2& x2,
const X3& x3,
const X4& x4,
const X5& x5,
double delta = 1e-5) {
477 "Template argument Y must be a manifold type.");
479 "Template argument X1 must be a manifold type.");
480 return numericalDerivative11<Y, X3>(boost::bind(h, boost::cref(x1), boost::cref(x2), _1, boost::cref(x4), boost::cref(x5)), x3, delta);
483 template<
class Y,
class X1,
class X2,
class X3,
class X4,
class X5>
484 inline typename internal::FixedSizeMatrix<Y,X3>::type
numericalDerivative53(Y (*h)(
const X1&,
const X2&,
const X3&,
const X4&,
const X5&),
485 const X1& x1,
const X2& x2,
const X3& x3,
const X4& x4,
const X5& x5,
double delta = 1e-5) {
486 return numericalDerivative53<Y, X1, X2, X3, X4, X5>(boost::bind(h, _1, _2, _3, _4, _5), x1, x2, x3, x4, x5);
500 template<
class Y,
class X1,
class X2,
class X3,
class X4,
class X5>
502 boost::function<Y(
const X1&,
const X2&,
const X3&,
const X4&,
const X5&)> h,
const X1& x1,
503 const X2& x2,
const X3& x3,
const X4& x4,
const X5& x5,
double delta = 1e-5) {
505 "Template argument Y must be a manifold type.");
507 "Template argument X1 must be a manifold type.");
508 return numericalDerivative11<Y, X4>(boost::bind(h, boost::cref(x1), boost::cref(x2), boost::cref(x3), _1, boost::cref(x5)), x4, delta);
511 template<
class Y,
class X1,
class X2,
class X3,
class X4,
class X5>
512 inline typename internal::FixedSizeMatrix<Y,X4>::type
numericalDerivative53(Y (*h)(
const X1&,
const X2&,
const X3&,
const X4&,
const X5&),
513 const X1& x1,
const X2& x2,
const X3& x3,
const X4& x4,
const X5& x5,
double delta = 1e-5) {
514 return numericalDerivative54<Y, X1, X2, X3, X4, X5>(boost::bind(h, _1, _2, _3, _4, _5), x1, x2, x3, x4, x5);
528 template<
class Y,
class X1,
class X2,
class X3,
class X4,
class X5>
530 boost::function<Y(
const X1&,
const X2&,
const X3&,
const X4&,
const X5&)> h,
const X1& x1,
531 const X2& x2,
const X3& x3,
const X4& x4,
const X5& x5,
double delta = 1e-5) {
533 "Template argument Y must be a manifold type.");
535 "Template argument X1 must be a manifold type.");
536 return numericalDerivative11<Y, X5>(boost::bind(h, boost::cref(x1), boost::cref(x2), boost::cref(x3), boost::cref(x4), _1), x5, delta);
539 template<
class Y,
class X1,
class X2,
class X3,
class X4,
class X5>
540 inline typename internal::FixedSizeMatrix<Y,X5>::type
numericalDerivative53(Y (*h)(
const X1&,
const X2&,
const X3&,
const X4&,
const X5&),
541 const X1& x1,
const X2& x2,
const X3& x3,
const X4& x4,
const X5& x5,
double delta = 1e-5) {
542 return numericalDerivative55<Y, X1, X2, X3, X4, X5>(boost::bind(h, _1, _2, _3, _4, _5), x1, x2, x3, x4, x5);
554 inline typename internal::FixedSizeMatrix<X,X>::type
numericalHessian(boost::function<
double(
const X&)> f,
const X& x,
555 double delta = 1e-5) {
557 "Template argument X must be a manifold type.");
558 typedef Eigen::Matrix<double, traits<X>::dimension, 1> VectorD;
559 typedef boost::function<double(
const X&)> F;
560 typedef boost::function<VectorD(F,
const X&,
double)> G;
561 G ng = static_cast<G>(numericalGradient<X> );
562 return numericalDerivative11<VectorD, X>(boost::bind(ng, f, _1, delta), x,
567 inline typename internal::FixedSizeMatrix<X,X>::type
numericalHessian(
double (*f)(
const X&),
const X& x,
double delta =
575 template<
class X1,
class X2>
577 const boost::function<double(
const X1&,
const X2&)>& f_;
581 typedef typename internal::FixedSizeMatrix<X1>::type Vector;
583 G_x1(
const boost::function<
double(
const X1&,
const X2&)>& f,
const X1& x1,
585 f_(f), x1_(x1), delta_(delta) {
587 Vector operator()(
const X2& x2) {
588 return numericalGradient<X1>(boost::bind(f_, _1, boost::cref(x2)), x1_, delta_);
592 template<
class X1,
class X2>
593 inline typename internal::FixedSizeMatrix<X1,X2>::type numericalHessian212(
594 boost::function<
double(
const X1&,
const X2&)> f,
const X1& x1,
const X2& x2,
595 double delta = 1e-5) {
596 typedef typename internal::FixedSizeMatrix<X1>::type Vector;
598 return numericalDerivative11<Vector, X2>(
599 boost::function<Vector(
const X2&)>(
600 boost::bind<Vector>(boost::ref(g_x1), _1)), x2, delta);
603 template<
class X1,
class X2>
604 inline typename internal::FixedSizeMatrix<X1,X2>::type numericalHessian212(
double (*f)(
const X1&,
const X2&),
605 const X1& x1,
const X2& x2,
double delta = 1e-5) {
606 return numericalHessian212(boost::function<
double(
const X1&,
const X2&)>(f),
610 template<
class X1,
class X2>
611 inline typename internal::FixedSizeMatrix<X1,X1>::type numericalHessian211(
612 boost::function<
double(
const X1&,
const X2&)> f,
const X1& x1,
const X2& x2,
613 double delta = 1e-5) {
615 typedef typename internal::FixedSizeMatrix<X1>::type Vector;
617 Vector (*numGrad)(boost::function<double(
const X1&)>,
const X1&,
618 double) = &numericalGradient<X1>;
619 boost::function<double(
const X1&)> f2(boost::bind(f, _1, boost::cref(x2)));
621 return numericalDerivative11<Vector, X1>(
622 boost::function<Vector(
const X1&)>(boost::bind(numGrad, f2, _1, delta)),
626 template<
class X1,
class X2>
627 inline typename internal::FixedSizeMatrix<X1,X1>::type numericalHessian211(
double (*f)(
const X1&,
const X2&),
628 const X1& x1,
const X2& x2,
double delta = 1e-5) {
629 return numericalHessian211(boost::function<
double(
const X1&,
const X2&)>(f),
633 template<
class X1,
class X2>
634 inline typename internal::FixedSizeMatrix<X2,X2>::type numericalHessian222(
635 boost::function<
double(
const X1&,
const X2&)> f,
const X1& x1,
const X2& x2,
636 double delta = 1e-5) {
637 typedef typename internal::FixedSizeMatrix<X2>::type Vector;
638 Vector (*numGrad)(boost::function<double(
const X2&)>,
const X2&,
639 double) = &numericalGradient<X2>;
640 boost::function<double(
const X2&)> f2(boost::bind(f, boost::cref(x1), _1));
642 return numericalDerivative11<Vector, X2>(
643 boost::function<Vector(
const X2&)>(boost::bind(numGrad, f2, _1, delta)),
647 template<
class X1,
class X2>
648 inline typename internal::FixedSizeMatrix<X2,X2>::type numericalHessian222(
double (*f)(
const X1&,
const X2&),
649 const X1& x1,
const X2& x2,
double delta = 1e-5) {
650 return numericalHessian222(boost::function<
double(
const X1&,
const X2&)>(f),
658 template<
class X1,
class X2,
class X3>
660 boost::function<
double(
const X1&,
const X2&,
const X3&)> f,
const X1& x1,
661 const X2& x2,
const X3& x3,
double delta = 1e-5) {
662 typedef typename internal::FixedSizeMatrix<X1>::type Vector;
663 Vector (*numGrad)(boost::function<double(
const X1&)>,
const X1&,
664 double) = &numericalGradient<X1>;
665 boost::function<double(
const X1&)> f2(boost::bind(f, _1, boost::cref(x2), boost::cref(x3)));
667 return numericalDerivative11<Vector, X1>(
668 boost::function<Vector(
const X1&)>(boost::bind(numGrad, f2, _1, delta)),
672 template<
class X1,
class X2,
class X3>
673 inline typename internal::FixedSizeMatrix<X1,X1>::type
numericalHessian311(
double (*f)(
const X1&,
const X2&,
const X3&),
674 const X1& x1,
const X2& x2,
const X3& x3,
double delta = 1e-5) {
676 boost::function<
double(
const X1&,
const X2&,
const X3&)>(f), x1, x2, x3,
681 template<
class X1,
class X2,
class X3>
682 inline typename internal::FixedSizeMatrix<X2,X2>::type numericalHessian322(
683 boost::function<
double(
const X1&,
const X2&,
const X3&)> f,
const X1& x1,
684 const X2& x2,
const X3& x3,
double delta = 1e-5) {
685 typedef typename internal::FixedSizeMatrix<X2>::type Vector;
686 Vector (*numGrad)(boost::function<double(
const X2&)>,
const X2&,
687 double) = &numericalGradient<X2>;
688 boost::function<double(
const X2&)> f2(boost::bind(f, boost::cref(x1), _1, boost::cref(x3)));
690 return numericalDerivative11<Vector, X2>(
691 boost::function<Vector(
const X2&)>(boost::bind(numGrad, f2, _1, delta)),
695 template<
class X1,
class X2,
class X3>
696 inline typename internal::FixedSizeMatrix<X2,X2>::type numericalHessian322(
double (*f)(
const X1&,
const X2&,
const X3&),
697 const X1& x1,
const X2& x2,
const X3& x3,
double delta = 1e-5) {
698 return numericalHessian322(
699 boost::function<
double(
const X1&,
const X2&,
const X3&)>(f), x1, x2, x3,
704 template<
class X1,
class X2,
class X3>
705 inline typename internal::FixedSizeMatrix<X3,X3>::type numericalHessian333(
706 boost::function<
double(
const X1&,
const X2&,
const X3&)> f,
const X1& x1,
707 const X2& x2,
const X3& x3,
double delta = 1e-5) {
708 typedef typename internal::FixedSizeMatrix<X3>::type Vector;
709 Vector (*numGrad)(boost::function<double(
const X3&)>,
const X3&,
710 double) = &numericalGradient<X3>;
711 boost::function<double(
const X3&)> f2(boost::bind(f, boost::cref(x1), boost::cref(x2), _1));
713 return numericalDerivative11<Vector, X3>(
714 boost::function<Vector(
const X3&)>(boost::bind(numGrad, f2, _1, delta)),
718 template<
class X1,
class X2,
class X3>
719 inline typename internal::FixedSizeMatrix<X3,X3>::type numericalHessian333(
double (*f)(
const X1&,
const X2&,
const X3&),
720 const X1& x1,
const X2& x2,
const X3& x3,
double delta = 1e-5) {
721 return numericalHessian333(
722 boost::function<
double(
const X1&,
const X2&,
const X3&)>(f), x1, x2, x3,
727 template<
class X1,
class X2,
class X3>
728 inline typename internal::FixedSizeMatrix<X1,X2>::type numericalHessian312(
729 boost::function<
double(
const X1&,
const X2&,
const X3&)> f,
const X1& x1,
730 const X2& x2,
const X3& x3,
double delta = 1e-5) {
731 return numericalHessian212<X1, X2>(
732 boost::function<
double(
const X1&,
const X2&)>(boost::bind(f, _1, _2, boost::cref(x3))),
736 template<
class X1,
class X2,
class X3>
737 inline typename internal::FixedSizeMatrix<X1,X3>::type numericalHessian313(
738 boost::function<
double(
const X1&,
const X2&,
const X3&)> f,
const X1& x1,
739 const X2& x2,
const X3& x3,
double delta = 1e-5) {
740 return numericalHessian212<X1, X3>(
741 boost::function<
double(
const X1&,
const X3&)>(boost::bind(f, _1, boost::cref(x2), _2)),
745 template<
class X1,
class X2,
class X3>
746 inline typename internal::FixedSizeMatrix<X2,X3>::type numericalHessian323(
747 boost::function<
double(
const X1&,
const X2&,
const X3&)> f,
const X1& x1,
748 const X2& x2,
const X3& x3,
double delta = 1e-5) {
749 return numericalHessian212<X2, X3>(
750 boost::function<
double(
const X2&,
const X3&)>(boost::bind(f, boost::cref(x1), _1, _2)),
755 template<
class X1,
class X2,
class X3>
756 inline typename internal::FixedSizeMatrix<X1,X2>::type numericalHessian312(
double (*f)(
const X1&,
const X2&,
const X3&),
757 const X1& x1,
const X2& x2,
const X3& x3,
double delta = 1e-5) {
758 return numericalHessian312(
759 boost::function<
double(
const X1&,
const X2&,
const X3&)>(f), x1, x2, x3,
763 template<
class X1,
class X2,
class X3>
764 inline typename internal::FixedSizeMatrix<X1,X3>::type numericalHessian313(
double (*f)(
const X1&,
const X2&,
const X3&),
765 const X1& x1,
const X2& x2,
const X3& x3,
double delta = 1e-5) {
766 return numericalHessian313(
767 boost::function<
double(
const X1&,
const X2&,
const X3&)>(f), x1, x2, x3,
771 template<
class X1,
class X2,
class X3>
772 inline typename internal::FixedSizeMatrix<X2,X3>::type numericalHessian323(
double (*f)(
const X1&,
const X2&,
const X3&),
773 const X1& x1,
const X2& x2,
const X3& x3,
double delta = 1e-5) {
774 return numericalHessian323(
775 boost::function<
double(
const X1&,
const X2&,
const X3&)>(f), x1, x2, x3,
internal::FixedSizeMatrix< Y, X1 >::type numericalDerivative51(boost::function< Y(const X1 &, const X2 &, const X3 &, const X4 &, const X5 &)> h, const X1 &x1, const X2 &x2, const X3 &x3, const X4 &x4, const X5 &x5, double delta=1e-5)
Compute numerical derivative in argument 1 of 5-argument function.
Definition: numericalDerivative.h:417
Definition: numericalDerivative.h:64
internal::FixedSizeMatrix< Y, X2 >::type numericalDerivative42(boost::function< Y(const X1 &, const X2 &, const X3 &, const X4 &)> h, const X1 &x1, const X2 &x2, const X3 &x3, const X4 &x4, double delta=1e-5)
Compute numerical derivative in argument 2 of 4-argument function.
Definition: numericalDerivative.h:335
internal::FixedSizeMatrix< Y, X1 >::type numericalDerivative31(boost::function< Y(const X1 &, const X2 &, const X3 &)> h, const X1 &x1, const X2 &x2, const X3 &x3, double delta=1e-5)
Compute numerical derivative in argument 1 of ternary function.
Definition: numericalDerivative.h:224
internal::FixedSizeMatrix< Y, X1 >::type numericalDerivative41(boost::function< Y(const X1 &, const X2 &, const X3 &, const X4 &)> h, const X1 &x1, const X2 &x2, const X3 &x3, const X4 &x4, double delta=1e-5)
Compute numerical derivative in argument 1 of 4-argument function.
Definition: numericalDerivative.h:308
internal::FixedSizeMatrix< Y, X1 >::type numericalDerivative21(const boost::function< Y(const X1 &, const X2 &)> &h, const X1 &x1, const X2 &x2, double delta=1e-5)
Compute numerical derivative in argument 1 of binary function.
Definition: numericalDerivative.h:172
Eigen::Matrix< double, N, 1 > numericalGradient(boost::function< double(const X &)> h, const X &x, double delta=1e-5)
Numerically compute gradient of scalar function.
Definition: numericalDerivative.h:78
internal::FixedSizeMatrix< Y, X4 >::type numericalDerivative54(boost::function< Y(const X1 &, const X2 &, const X3 &, const X4 &, const X5 &)> h, const X1 &x1, const X2 &x2, const X3 &x3, const X4 &x4, const X5 &x5, double delta=1e-5)
Compute numerical derivative in argument 4 of 5-argument function.
Definition: numericalDerivative.h:501
internal::FixedSizeMatrix< Y, X5 >::type numericalDerivative55(boost::function< Y(const X1 &, const X2 &, const X3 &, const X4 &, const X5 &)> h, const X1 &x1, const X2 &x2, const X3 &x3, const X4 &x4, const X5 &x5, double delta=1e-5)
Compute numerical derivative in argument 5 of 5-argument function.
Definition: numericalDerivative.h:529
internal::FixedSizeMatrix< X1, X1 >::type numericalHessian311(boost::function< double(const X1 &, const X2 &, const X3 &)> f, const X1 &x1, const X2 &x2, const X3 &x3, double delta=1e-5)
Numerical Hessian for tenary functions.
Definition: numericalDerivative.h:659
internal::FixedSizeMatrix< X, X >::type numericalHessian(boost::function< double(const X &)> f, const X &x, double delta=1e-5)
Compute numerical Hessian matrix.
Definition: numericalDerivative.h:554
Helper class that computes the derivative of f w.r.t.
Definition: numericalDerivative.h:576
internal::FixedSizeMatrix< Y, X2 >::type numericalDerivative22(boost::function< Y(const X1 &, const X2 &)> h, const X1 &x1, const X2 &x2, double delta=1e-5)
Compute numerical derivative in argument 2 of binary function.
Definition: numericalDerivative.h:197
tag to assert a type is a manifold
Definition: Manifold.h:33
A manifold defines a space in which there is a notion of a linear tangent space that can be centered ...
Definition: concepts.h:30
A non-templated config holding any types of Manifold-group elements.
internal::FixedSizeMatrix< Y, X3 >::type numericalDerivative43(boost::function< Y(const X1 &, const X2 &, const X3 &, const X4 &)> h, const X1 &x1, const X2 &x2, const X3 &x3, const X4 &x4, double delta=1e-5)
Compute numerical derivative in argument 3 of 4-argument function.
Definition: numericalDerivative.h:362
internal::FixedSizeMatrix< Y, X2 >::type numericalDerivative32(boost::function< Y(const X1 &, const X2 &, const X3 &)> h, const X1 &x1, const X2 &x2, const X3 &x3, double delta=1e-5)
Compute numerical derivative in argument 2 of ternary function.
Definition: numericalDerivative.h:252
Global functions in a separate testing namespace.
Definition: chartTesting.h:28
internal::FixedSizeMatrix< Y, X >::type numericalDerivative11(boost::function< Y(const X &)> h, const X &x, double delta=1e-5)
New-style numerical derivatives using manifold_traits.
Definition: numericalDerivative.h:118
internal::FixedSizeMatrix< Y, X2 >::type numericalDerivative52(boost::function< Y(const X1 &, const X2 &, const X3 &, const X4 &, const X5 &)> h, const X1 &x1, const X2 &x2, const X3 &x3, const X4 &x4, const X5 &x5, double delta=1e-5)
Compute numerical derivative in argument 2 of 5-argument function.
Definition: numericalDerivative.h:445
internal::FixedSizeMatrix< Y, X3 >::type numericalDerivative53(boost::function< Y(const X1 &, const X2 &, const X3 &, const X4 &, const X5 &)> h, const X1 &x1, const X2 &x2, const X3 &x3, const X4 &x4, const X5 &x5, double delta=1e-5)
Compute numerical derivative in argument 3 of 5-argument function.
Definition: numericalDerivative.h:473
internal::FixedSizeMatrix< Y, X4 >::type numericalDerivative44(boost::function< Y(const X1 &, const X2 &, const X3 &, const X4 &)> h, const X1 &x1, const X2 &x2, const X3 &x3, const X4 &x4, double delta=1e-5)
Compute numerical derivative in argument 4 of 4-argument function.
Definition: numericalDerivative.h:389
Base class and basic functions for Lie types.
internal::FixedSizeMatrix< Y, X3 >::type numericalDerivative33(boost::function< Y(const X1 &, const X2 &, const X3 &)> h, const X1 &x1, const X2 &x2, const X3 &x3, double delta=1e-5)
Compute numerical derivative in argument 3 of ternary function.
Definition: numericalDerivative.h:280