Imported Upstream version ceres 1.13.0
[platform/upstream/ceres-solver.git] / include / ceres / jet.h
1 // Ceres Solver - A fast non-linear least squares minimizer
2 // Copyright 2015 Google Inc. All rights reserved.
3 // http://ceres-solver.org/
4 //
5 // Redistribution and use in source and binary forms, with or without
6 // modification, are permitted provided that the following conditions are met:
7 //
8 // * Redistributions of source code must retain the above copyright notice,
9 //   this list of conditions and the following disclaimer.
10 // * Redistributions in binary form must reproduce the above copyright notice,
11 //   this list of conditions and the following disclaimer in the documentation
12 //   and/or other materials provided with the distribution.
13 // * Neither the name of Google Inc. nor the names of its contributors may be
14 //   used to endorse or promote products derived from this software without
15 //   specific prior written permission.
16 //
17 // THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
18 // AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
19 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
20 // ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE
21 // LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
22 // CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
23 // SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
24 // INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
25 // CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
26 // ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
27 // POSSIBILITY OF SUCH DAMAGE.
28 //
29 // Author: keir@google.com (Keir Mierle)
30 //
31 // A simple implementation of N-dimensional dual numbers, for automatically
32 // computing exact derivatives of functions.
33 //
34 // While a complete treatment of the mechanics of automatic differentation is
35 // beyond the scope of this header (see
36 // http://en.wikipedia.org/wiki/Automatic_differentiation for details), the
37 // basic idea is to extend normal arithmetic with an extra element, "e," often
38 // denoted with the greek symbol epsilon, such that e != 0 but e^2 = 0. Dual
39 // numbers are extensions of the real numbers analogous to complex numbers:
40 // whereas complex numbers augment the reals by introducing an imaginary unit i
41 // such that i^2 = -1, dual numbers introduce an "infinitesimal" unit e such
42 // that e^2 = 0. Dual numbers have two components: the "real" component and the
43 // "infinitesimal" component, generally written as x + y*e. Surprisingly, this
44 // leads to a convenient method for computing exact derivatives without needing
45 // to manipulate complicated symbolic expressions.
46 //
47 // For example, consider the function
48 //
49 //   f(x) = x^2 ,
50 //
51 // evaluated at 10. Using normal arithmetic, f(10) = 100, and df/dx(10) = 20.
52 // Next, augument 10 with an infinitesimal to get:
53 //
54 //   f(10 + e) = (10 + e)^2
55 //             = 100 + 2 * 10 * e + e^2
56 //             = 100 + 20 * e       -+-
57 //                     --            |
58 //                     |             +--- This is zero, since e^2 = 0
59 //                     |
60 //                     +----------------- This is df/dx!
61 //
62 // Note that the derivative of f with respect to x is simply the infinitesimal
63 // component of the value of f(x + e). So, in order to take the derivative of
64 // any function, it is only necessary to replace the numeric "object" used in
65 // the function with one extended with infinitesimals. The class Jet, defined in
66 // this header, is one such example of this, where substitution is done with
67 // templates.
68 //
69 // To handle derivatives of functions taking multiple arguments, different
70 // infinitesimals are used, one for each variable to take the derivative of. For
71 // example, consider a scalar function of two scalar parameters x and y:
72 //
73 //   f(x, y) = x^2 + x * y
74 //
75 // Following the technique above, to compute the derivatives df/dx and df/dy for
76 // f(1, 3) involves doing two evaluations of f, the first time replacing x with
77 // x + e, the second time replacing y with y + e.
78 //
79 // For df/dx:
80 //
81 //   f(1 + e, y) = (1 + e)^2 + (1 + e) * 3
82 //               = 1 + 2 * e + 3 + 3 * e
83 //               = 4 + 5 * e
84 //
85 //               --> df/dx = 5
86 //
87 // For df/dy:
88 //
89 //   f(1, 3 + e) = 1^2 + 1 * (3 + e)
90 //               = 1 + 3 + e
91 //               = 4 + e
92 //
93 //               --> df/dy = 1
94 //
95 // To take the gradient of f with the implementation of dual numbers ("jets") in
96 // this file, it is necessary to create a single jet type which has components
97 // for the derivative in x and y, and passing them to a templated version of f:
98 //
99 //   template<typename T>
100 //   T f(const T &x, const T &y) {
101 //     return x * x + x * y;
102 //   }
103 //
104 //   // The "2" means there should be 2 dual number components.
105 //   Jet<double, 2> x(0);  // Pick the 0th dual number for x.
106 //   Jet<double, 2> y(1);  // Pick the 1st dual number for y.
107 //   Jet<double, 2> z = f(x, y);
108 //
109 //   LOG(INFO) << "df/dx = " << z.v[0]
110 //             << "df/dy = " << z.v[1];
111 //
112 // Most users should not use Jet objects directly; a wrapper around Jet objects,
113 // which makes computing the derivative, gradient, or jacobian of templated
114 // functors simple, is in autodiff.h. Even autodiff.h should not be used
115 // directly; instead autodiff_cost_function.h is typically the file of interest.
116 //
117 // For the more mathematically inclined, this file implements first-order
118 // "jets". A 1st order jet is an element of the ring
119 //
120 //   T[N] = T[t_1, ..., t_N] / (t_1, ..., t_N)^2
121 //
122 // which essentially means that each jet consists of a "scalar" value 'a' from T
123 // and a 1st order perturbation vector 'v' of length N:
124 //
125 //   x = a + \sum_i v[i] t_i
126 //
127 // A shorthand is to write an element as x = a + u, where u is the pertubation.
128 // Then, the main point about the arithmetic of jets is that the product of
129 // perturbations is zero:
130 //
131 //   (a + u) * (b + v) = ab + av + bu + uv
132 //                     = ab + (av + bu) + 0
133 //
134 // which is what operator* implements below. Addition is simpler:
135 //
136 //   (a + u) + (b + v) = (a + b) + (u + v).
137 //
138 // The only remaining question is how to evaluate the function of a jet, for
139 // which we use the chain rule:
140 //
141 //   f(a + u) = f(a) + f'(a) u
142 //
143 // where f'(a) is the (scalar) derivative of f at a.
144 //
145 // By pushing these things through sufficiently and suitably templated
146 // functions, we can do automatic differentiation. Just be sure to turn on
147 // function inlining and common-subexpression elimination, or it will be very
148 // slow!
149 //
150 // WARNING: Most Ceres users should not directly include this file or know the
151 // details of how jets work. Instead the suggested method for automatic
152 // derivatives is to use autodiff_cost_function.h, which is a wrapper around
153 // both jets.h and autodiff.h to make taking derivatives of cost functions for
154 // use in Ceres easier.
155
156 #ifndef CERES_PUBLIC_JET_H_
157 #define CERES_PUBLIC_JET_H_
158
159 #include <cmath>
160 #include <iosfwd>
161 #include <iostream>  // NOLINT
162 #include <limits>
163 #include <string>
164
165 #include "Eigen/Core"
166 #include "ceres/fpclassify.h"
167 #include "ceres/internal/port.h"
168
169 namespace ceres {
170
171 template <typename T, int N>
172 struct Jet {
173   enum { DIMENSION = N };
174
175   // Default-construct "a" because otherwise this can lead to false errors about
176   // uninitialized uses when other classes relying on default constructed T
177   // (where T is a Jet<T, N>). This usually only happens in opt mode. Note that
178   // the C++ standard mandates that e.g. default constructed doubles are
179   // initialized to 0.0; see sections 8.5 of the C++03 standard.
180   Jet() : a() {
181     v.setZero();
182   }
183
184   // Constructor from scalar: a + 0.
185   explicit Jet(const T& value) {
186     a = value;
187     v.setZero();
188   }
189
190   // Constructor from scalar plus variable: a + t_i.
191   Jet(const T& value, int k) {
192     a = value;
193     v.setZero();
194     v[k] = T(1.0);
195   }
196
197   // Constructor from scalar and vector part
198   // The use of Eigen::DenseBase allows Eigen expressions
199   // to be passed in without being fully evaluated until
200   // they are assigned to v
201   template<typename Derived>
202   EIGEN_STRONG_INLINE Jet(const T& a, const Eigen::DenseBase<Derived> &v)
203       : a(a), v(v) {
204   }
205
206   // Compound operators
207   Jet<T, N>& operator+=(const Jet<T, N> &y) {
208     *this = *this + y;
209     return *this;
210   }
211
212   Jet<T, N>& operator-=(const Jet<T, N> &y) {
213     *this = *this - y;
214     return *this;
215   }
216
217   Jet<T, N>& operator*=(const Jet<T, N> &y) {
218     *this = *this * y;
219     return *this;
220   }
221
222   Jet<T, N>& operator/=(const Jet<T, N> &y) {
223     *this = *this / y;
224     return *this;
225   }
226
227   // Compound with scalar operators.
228   Jet<T, N>& operator+=(const T& s) {
229     *this = *this + s;
230     return *this;
231   }
232
233   Jet<T, N>& operator-=(const T& s) {
234     *this = *this - s;
235     return *this;
236   }
237
238   Jet<T, N>& operator*=(const T& s) {
239     *this = *this * s;
240     return *this;
241   }
242
243   Jet<T, N>& operator/=(const T& s) {
244     *this = *this / s;
245     return *this;
246   }
247
248   // The scalar part.
249   T a;
250
251   // The infinitesimal part.
252   //
253   // We allocate Jets on the stack and other places they might not be aligned
254   // to X(=16 [SSE], 32 [AVX] etc)-byte boundaries, which would prevent the safe
255   // use of vectorisation.  If we have C++11, we can specify the alignment.
256   // However, the standard gives wide lattitude as to what alignments are valid,
257   // and it might be that the maximum supported alignment *guaranteed* to be
258   // supported is < 16, in which case we do not specify an alignment, as this
259   // implies the host is not a modern x86 machine.  If using < C++11, we cannot
260   // specify alignment.
261 #ifndef CERES_USE_CXX11
262   // Without >= C++11, we cannot specify the alignment so fall back to safe,
263   // unvectorised version.
264   Eigen::Matrix<T, N, 1, Eigen::DontAlign> v;
265 #else
266   // Enable vectorisation iff the maximum supported scalar alignment is >=
267   // 16 bytes, as this is the minimum required by Eigen for any vectorisation.
268   //
269   // NOTE: It might be the case that we could get >= 16-byte alignment even if
270   //       kMaxAlignBytes < 16.  However we can't guarantee that this
271   //       would happen (and it should not for any modern x86 machine) and if it
272   //       didn't, we could get misaligned Jets.
273   static constexpr int kAlignOrNot =
274       16 <= ::ceres::port_constants::kMaxAlignBytes
275             ? Eigen::AutoAlign : Eigen::DontAlign;
276 #if defined(EIGEN_MAX_ALIGN_BYTES)
277   // Eigen >= 3.3 supports AVX & FMA instructions that require 32-byte alignment
278   // (greater for AVX512).  Rather than duplicating the detection logic, use
279   // Eigen's macro for the alignment size.
280   //
281   // NOTE: EIGEN_MAX_ALIGN_BYTES can be > 16 (e.g. 32 for AVX), even though
282   //       kMaxAlignBytes will max out at 16.  We are therefore relying on
283   //       Eigen's detection logic to ensure that this does not result in
284   //       misaligned Jets.
285 #define CERES_JET_ALIGN_BYTES EIGEN_MAX_ALIGN_BYTES
286 #else
287   // Eigen < 3.3 only supported 16-byte alignment.
288 #define CERES_JET_ALIGN_BYTES 16
289 #endif
290   // Default to the native alignment if 16-byte alignment is not guaranteed to
291   // be supported.  We cannot use alignof(T) as if we do, GCC 4.8 complains that
292   // the alignment 'is not an integer constant', although Clang accepts it.
293   static constexpr size_t kAlignment = kAlignOrNot == Eigen::AutoAlign
294             ? CERES_JET_ALIGN_BYTES : alignof(double);
295 #undef CERES_JET_ALIGN_BYTES
296   alignas(kAlignment) Eigen::Matrix<T, N, 1, kAlignOrNot> v;
297 #endif
298 };
299
300 // Unary +
301 template<typename T, int N> inline
302 Jet<T, N> const& operator+(const Jet<T, N>& f) {
303   return f;
304 }
305
306 // TODO(keir): Try adding __attribute__((always_inline)) to these functions to
307 // see if it causes a performance increase.
308
309 // Unary -
310 template<typename T, int N> inline
311 Jet<T, N> operator-(const Jet<T, N>&f) {
312   return Jet<T, N>(-f.a, -f.v);
313 }
314
315 // Binary +
316 template<typename T, int N> inline
317 Jet<T, N> operator+(const Jet<T, N>& f,
318                     const Jet<T, N>& g) {
319   return Jet<T, N>(f.a + g.a, f.v + g.v);
320 }
321
322 // Binary + with a scalar: x + s
323 template<typename T, int N> inline
324 Jet<T, N> operator+(const Jet<T, N>& f, T s) {
325   return Jet<T, N>(f.a + s, f.v);
326 }
327
328 // Binary + with a scalar: s + x
329 template<typename T, int N> inline
330 Jet<T, N> operator+(T s, const Jet<T, N>& f) {
331   return Jet<T, N>(f.a + s, f.v);
332 }
333
334 // Binary -
335 template<typename T, int N> inline
336 Jet<T, N> operator-(const Jet<T, N>& f,
337                     const Jet<T, N>& g) {
338   return Jet<T, N>(f.a - g.a, f.v - g.v);
339 }
340
341 // Binary - with a scalar: x - s
342 template<typename T, int N> inline
343 Jet<T, N> operator-(const Jet<T, N>& f, T s) {
344   return Jet<T, N>(f.a - s, f.v);
345 }
346
347 // Binary - with a scalar: s - x
348 template<typename T, int N> inline
349 Jet<T, N> operator-(T s, const Jet<T, N>& f) {
350   return Jet<T, N>(s - f.a, -f.v);
351 }
352
353 // Binary *
354 template<typename T, int N> inline
355 Jet<T, N> operator*(const Jet<T, N>& f,
356                     const Jet<T, N>& g) {
357   return Jet<T, N>(f.a * g.a, f.a * g.v + f.v * g.a);
358 }
359
360 // Binary * with a scalar: x * s
361 template<typename T, int N> inline
362 Jet<T, N> operator*(const Jet<T, N>& f, T s) {
363   return Jet<T, N>(f.a * s, f.v * s);
364 }
365
366 // Binary * with a scalar: s * x
367 template<typename T, int N> inline
368 Jet<T, N> operator*(T s, const Jet<T, N>& f) {
369   return Jet<T, N>(f.a * s, f.v * s);
370 }
371
372 // Binary /
373 template<typename T, int N> inline
374 Jet<T, N> operator/(const Jet<T, N>& f,
375                     const Jet<T, N>& g) {
376   // This uses:
377   //
378   //   a + u   (a + u)(b - v)   (a + u)(b - v)
379   //   ----- = -------------- = --------------
380   //   b + v   (b + v)(b - v)        b^2
381   //
382   // which holds because v*v = 0.
383   const T g_a_inverse = T(1.0) / g.a;
384   const T f_a_by_g_a = f.a * g_a_inverse;
385   return Jet<T, N>(f.a * g_a_inverse, (f.v - f_a_by_g_a * g.v) * g_a_inverse);
386 }
387
388 // Binary / with a scalar: s / x
389 template<typename T, int N> inline
390 Jet<T, N> operator/(T s, const Jet<T, N>& g) {
391   const T minus_s_g_a_inverse2 = -s / (g.a * g.a);
392   return Jet<T, N>(s / g.a, g.v * minus_s_g_a_inverse2);
393 }
394
395 // Binary / with a scalar: x / s
396 template<typename T, int N> inline
397 Jet<T, N> operator/(const Jet<T, N>& f, T s) {
398   const T s_inverse = T(1.0) / s;
399   return Jet<T, N>(f.a * s_inverse, f.v * s_inverse);
400 }
401
402 // Binary comparison operators for both scalars and jets.
403 #define CERES_DEFINE_JET_COMPARISON_OPERATOR(op) \
404 template<typename T, int N> inline \
405 bool operator op(const Jet<T, N>& f, const Jet<T, N>& g) { \
406   return f.a op g.a; \
407 } \
408 template<typename T, int N> inline \
409 bool operator op(const T& s, const Jet<T, N>& g) { \
410   return s op g.a; \
411 } \
412 template<typename T, int N> inline \
413 bool operator op(const Jet<T, N>& f, const T& s) { \
414   return f.a op s; \
415 }
416 CERES_DEFINE_JET_COMPARISON_OPERATOR( <  )  // NOLINT
417 CERES_DEFINE_JET_COMPARISON_OPERATOR( <= )  // NOLINT
418 CERES_DEFINE_JET_COMPARISON_OPERATOR( >  )  // NOLINT
419 CERES_DEFINE_JET_COMPARISON_OPERATOR( >= )  // NOLINT
420 CERES_DEFINE_JET_COMPARISON_OPERATOR( == )  // NOLINT
421 CERES_DEFINE_JET_COMPARISON_OPERATOR( != )  // NOLINT
422 #undef CERES_DEFINE_JET_COMPARISON_OPERATOR
423
424 // Pull some functions from namespace std.
425 //
426 // This is necessary because we want to use the same name (e.g. 'sqrt') for
427 // double-valued and Jet-valued functions, but we are not allowed to put
428 // Jet-valued functions inside namespace std.
429 //
430 // TODO(keir): Switch to "using".
431 inline double abs     (double x) { return std::abs(x);      }
432 inline double log     (double x) { return std::log(x);      }
433 inline double exp     (double x) { return std::exp(x);      }
434 inline double sqrt    (double x) { return std::sqrt(x);     }
435 inline double cos     (double x) { return std::cos(x);      }
436 inline double acos    (double x) { return std::acos(x);     }
437 inline double sin     (double x) { return std::sin(x);      }
438 inline double asin    (double x) { return std::asin(x);     }
439 inline double tan     (double x) { return std::tan(x);      }
440 inline double atan    (double x) { return std::atan(x);     }
441 inline double sinh    (double x) { return std::sinh(x);     }
442 inline double cosh    (double x) { return std::cosh(x);     }
443 inline double tanh    (double x) { return std::tanh(x);     }
444 inline double floor   (double x) { return std::floor(x);    }
445 inline double ceil    (double x) { return std::ceil(x);     }
446 inline double pow  (double x, double y) { return std::pow(x, y);   }
447 inline double atan2(double y, double x) { return std::atan2(y, x); }
448
449 // In general, f(a + h) ~= f(a) + f'(a) h, via the chain rule.
450
451 // abs(x + h) ~= x + h or -(x + h)
452 template <typename T, int N> inline
453 Jet<T, N> abs(const Jet<T, N>& f) {
454   return f.a < T(0.0) ? -f : f;
455 }
456
457 // log(a + h) ~= log(a) + h / a
458 template <typename T, int N> inline
459 Jet<T, N> log(const Jet<T, N>& f) {
460   const T a_inverse = T(1.0) / f.a;
461   return Jet<T, N>(log(f.a), f.v * a_inverse);
462 }
463
464 // exp(a + h) ~= exp(a) + exp(a) h
465 template <typename T, int N> inline
466 Jet<T, N> exp(const Jet<T, N>& f) {
467   const T tmp = exp(f.a);
468   return Jet<T, N>(tmp, tmp * f.v);
469 }
470
471 // sqrt(a + h) ~= sqrt(a) + h / (2 sqrt(a))
472 template <typename T, int N> inline
473 Jet<T, N> sqrt(const Jet<T, N>& f) {
474   const T tmp = sqrt(f.a);
475   const T two_a_inverse = T(1.0) / (T(2.0) * tmp);
476   return Jet<T, N>(tmp, f.v * two_a_inverse);
477 }
478
479 // cos(a + h) ~= cos(a) - sin(a) h
480 template <typename T, int N> inline
481 Jet<T, N> cos(const Jet<T, N>& f) {
482   return Jet<T, N>(cos(f.a), - sin(f.a) * f.v);
483 }
484
485 // acos(a + h) ~= acos(a) - 1 / sqrt(1 - a^2) h
486 template <typename T, int N> inline
487 Jet<T, N> acos(const Jet<T, N>& f) {
488   const T tmp = - T(1.0) / sqrt(T(1.0) - f.a * f.a);
489   return Jet<T, N>(acos(f.a), tmp * f.v);
490 }
491
492 // sin(a + h) ~= sin(a) + cos(a) h
493 template <typename T, int N> inline
494 Jet<T, N> sin(const Jet<T, N>& f) {
495   return Jet<T, N>(sin(f.a), cos(f.a) * f.v);
496 }
497
498 // asin(a + h) ~= asin(a) + 1 / sqrt(1 - a^2) h
499 template <typename T, int N> inline
500 Jet<T, N> asin(const Jet<T, N>& f) {
501   const T tmp = T(1.0) / sqrt(T(1.0) - f.a * f.a);
502   return Jet<T, N>(asin(f.a), tmp * f.v);
503 }
504
505 // tan(a + h) ~= tan(a) + (1 + tan(a)^2) h
506 template <typename T, int N> inline
507 Jet<T, N> tan(const Jet<T, N>& f) {
508   const T tan_a = tan(f.a);
509   const T tmp = T(1.0) + tan_a * tan_a;
510   return Jet<T, N>(tan_a, tmp * f.v);
511 }
512
513 // atan(a + h) ~= atan(a) + 1 / (1 + a^2) h
514 template <typename T, int N> inline
515 Jet<T, N> atan(const Jet<T, N>& f) {
516   const T tmp = T(1.0) / (T(1.0) + f.a * f.a);
517   return Jet<T, N>(atan(f.a), tmp * f.v);
518 }
519
520 // sinh(a + h) ~= sinh(a) + cosh(a) h
521 template <typename T, int N> inline
522 Jet<T, N> sinh(const Jet<T, N>& f) {
523   return Jet<T, N>(sinh(f.a), cosh(f.a) * f.v);
524 }
525
526 // cosh(a + h) ~= cosh(a) + sinh(a) h
527 template <typename T, int N> inline
528 Jet<T, N> cosh(const Jet<T, N>& f) {
529   return Jet<T, N>(cosh(f.a), sinh(f.a) * f.v);
530 }
531
532 // tanh(a + h) ~= tanh(a) + (1 - tanh(a)^2) h
533 template <typename T, int N> inline
534 Jet<T, N> tanh(const Jet<T, N>& f) {
535   const T tanh_a = tanh(f.a);
536   const T tmp = T(1.0) - tanh_a * tanh_a;
537   return Jet<T, N>(tanh_a, tmp * f.v);
538 }
539
540 // The floor function should be used with extreme care as this operation will
541 // result in a zero derivative which provides no information to the solver.
542 //
543 // floor(a + h) ~= floor(a) + 0
544 template <typename T, int N> inline
545 Jet<T, N> floor(const Jet<T, N>& f) {
546   return Jet<T, N>(floor(f.a));
547 }
548
549 // The ceil function should be used with extreme care as this operation will
550 // result in a zero derivative which provides no information to the solver.
551 //
552 // ceil(a + h) ~= ceil(a) + 0
553 template <typename T, int N> inline
554 Jet<T, N> ceil(const Jet<T, N>& f) {
555   return Jet<T, N>(ceil(f.a));
556 }
557
558 // Bessel functions of the first kind with integer order equal to 0, 1, n.
559 //
560 // Microsoft has deprecated the j[0,1,n]() POSIX Bessel functions in favour of
561 // _j[0,1,n]().  Where available on MSVC, use _j[0,1,n]() to avoid deprecated
562 // function errors in client code (the specific warning is suppressed when
563 // Ceres itself is built).
564 inline double BesselJ0(double x) {
565 #if defined(CERES_MSVC_USE_UNDERSCORE_PREFIXED_BESSEL_FUNCTIONS)
566   return _j0(x);
567 #else
568   return j0(x);
569 #endif
570 }
571 inline double BesselJ1(double x) {
572 #if defined(CERES_MSVC_USE_UNDERSCORE_PREFIXED_BESSEL_FUNCTIONS)
573   return _j1(x);
574 #else
575   return j1(x);
576 #endif
577 }
578 inline double BesselJn(int n, double x) {
579 #if defined(CERES_MSVC_USE_UNDERSCORE_PREFIXED_BESSEL_FUNCTIONS)
580   return _jn(n, x);
581 #else
582   return jn(n, x);
583 #endif
584 }
585
586 // For the formulae of the derivatives of the Bessel functions see the book:
587 // Olver, Lozier, Boisvert, Clark, NIST Handbook of Mathematical Functions,
588 // Cambridge University Press 2010.
589 //
590 // Formulae are also available at http://dlmf.nist.gov
591
592 // See formula http://dlmf.nist.gov/10.6#E3
593 // j0(a + h) ~= j0(a) - j1(a) h
594 template <typename T, int N> inline
595 Jet<T, N> BesselJ0(const Jet<T, N>& f) {
596   return Jet<T, N>(BesselJ0(f.a),
597                    -BesselJ1(f.a) * f.v);
598 }
599
600 // See formula http://dlmf.nist.gov/10.6#E1
601 // j1(a + h) ~= j1(a) + 0.5 ( j0(a) - j2(a) ) h
602 template <typename T, int N> inline
603 Jet<T, N> BesselJ1(const Jet<T, N>& f) {
604   return Jet<T, N>(BesselJ1(f.a),
605                    T(0.5) * (BesselJ0(f.a) - BesselJn(2, f.a)) * f.v);
606 }
607
608 // See formula http://dlmf.nist.gov/10.6#E1
609 // j_n(a + h) ~= j_n(a) + 0.5 ( j_{n-1}(a) - j_{n+1}(a) ) h
610 template <typename T, int N> inline
611 Jet<T, N> BesselJn(int n, const Jet<T, N>& f) {
612   return Jet<T, N>(BesselJn(n, f.a),
613                    T(0.5) * (BesselJn(n - 1, f.a) - BesselJn(n + 1, f.a)) * f.v);
614 }
615
616 // Jet Classification. It is not clear what the appropriate semantics are for
617 // these classifications. This picks that IsFinite and isnormal are "all"
618 // operations, i.e. all elements of the jet must be finite for the jet itself
619 // to be finite (or normal). For IsNaN and IsInfinite, the answer is less
620 // clear. This takes a "any" approach for IsNaN and IsInfinite such that if any
621 // part of a jet is nan or inf, then the entire jet is nan or inf. This leads
622 // to strange situations like a jet can be both IsInfinite and IsNaN, but in
623 // practice the "any" semantics are the most useful for e.g. checking that
624 // derivatives are sane.
625
626 // The jet is finite if all parts of the jet are finite.
627 template <typename T, int N> inline
628 bool IsFinite(const Jet<T, N>& f) {
629   if (!IsFinite(f.a)) {
630     return false;
631   }
632   for (int i = 0; i < N; ++i) {
633     if (!IsFinite(f.v[i])) {
634       return false;
635     }
636   }
637   return true;
638 }
639
640 // The jet is infinite if any part of the jet is infinite.
641 template <typename T, int N> inline
642 bool IsInfinite(const Jet<T, N>& f) {
643   if (IsInfinite(f.a)) {
644     return true;
645   }
646   for (int i = 0; i < N; i++) {
647     if (IsInfinite(f.v[i])) {
648       return true;
649     }
650   }
651   return false;
652 }
653
654 // The jet is NaN if any part of the jet is NaN.
655 template <typename T, int N> inline
656 bool IsNaN(const Jet<T, N>& f) {
657   if (IsNaN(f.a)) {
658     return true;
659   }
660   for (int i = 0; i < N; ++i) {
661     if (IsNaN(f.v[i])) {
662       return true;
663     }
664   }
665   return false;
666 }
667
668 // The jet is normal if all parts of the jet are normal.
669 template <typename T, int N> inline
670 bool IsNormal(const Jet<T, N>& f) {
671   if (!IsNormal(f.a)) {
672     return false;
673   }
674   for (int i = 0; i < N; ++i) {
675     if (!IsNormal(f.v[i])) {
676       return false;
677     }
678   }
679   return true;
680 }
681
682 // atan2(b + db, a + da) ~= atan2(b, a) + (- b da + a db) / (a^2 + b^2)
683 //
684 // In words: the rate of change of theta is 1/r times the rate of
685 // change of (x, y) in the positive angular direction.
686 template <typename T, int N> inline
687 Jet<T, N> atan2(const Jet<T, N>& g, const Jet<T, N>& f) {
688   // Note order of arguments:
689   //
690   //   f = a + da
691   //   g = b + db
692
693   T const tmp = T(1.0) / (f.a * f.a + g.a * g.a);
694   return Jet<T, N>(atan2(g.a, f.a), tmp * (- g.a * f.v + f.a * g.v));
695 }
696
697
698 // pow -- base is a differentiable function, exponent is a constant.
699 // (a+da)^p ~= a^p + p*a^(p-1) da
700 template <typename T, int N> inline
701 Jet<T, N> pow(const Jet<T, N>& f, double g) {
702   T const tmp = g * pow(f.a, g - T(1.0));
703   return Jet<T, N>(pow(f.a, g), tmp * f.v);
704 }
705
706 // pow -- base is a constant, exponent is a differentiable function.
707 // We have various special cases, see the comment for pow(Jet, Jet) for
708 // analysis:
709 //
710 // 1. For f > 0 we have: (f)^(g + dg) ~= f^g + f^g log(f) dg
711 //
712 // 2. For f == 0 and g > 0 we have: (f)^(g + dg) ~= f^g
713 //
714 // 3. For f < 0 and integer g we have: (f)^(g + dg) ~= f^g but if dg
715 // != 0, the derivatives are not defined and we return NaN.
716
717 template <typename T, int N> inline
718 Jet<T, N> pow(double f, const Jet<T, N>& g) {
719   if (f == 0 && g.a > 0) {
720     // Handle case 2.
721     return Jet<T, N>(T(0.0));
722   }
723   if (f < 0 && g.a == floor(g.a)) {
724     // Handle case 3.
725     Jet<T, N> ret(pow(f, g.a));
726     for (int i = 0; i < N; i++) {
727       if (g.v[i] != T(0.0)) {
728         // Return a NaN when g.v != 0.
729         ret.v[i] = std::numeric_limits<T>::quiet_NaN();
730       }
731     }
732     return ret;
733   }
734   // Handle case 1.
735   T const tmp = pow(f, g.a);
736   return Jet<T, N>(tmp, log(f) * tmp * g.v);
737 }
738
739 // pow -- both base and exponent are differentiable functions. This has a
740 // variety of special cases that require careful handling.
741 //
742 // 1. For f > 0:
743 //    (f + df)^(g + dg) ~= f^g + f^(g - 1) * (g * df + f * log(f) * dg)
744 //    The numerical evaluation of f * log(f) for f > 0 is well behaved, even for
745 //    extremely small values (e.g. 1e-99).
746 //
747 // 2. For f == 0 and g > 1: (f + df)^(g + dg) ~= 0
748 //    This cases is needed because log(0) can not be evaluated in the f > 0
749 //    expression. However the function f*log(f) is well behaved around f == 0
750 //    and its limit as f-->0 is zero.
751 //
752 // 3. For f == 0 and g == 1: (f + df)^(g + dg) ~= 0 + df
753 //
754 // 4. For f == 0 and 0 < g < 1: The value is finite but the derivatives are not.
755 //
756 // 5. For f == 0 and g < 0: The value and derivatives of f^g are not finite.
757 //
758 // 6. For f == 0 and g == 0: The C standard incorrectly defines 0^0 to be 1
759 //    "because there are applications that can exploit this definition". We
760 //    (arbitrarily) decree that derivatives here will be nonfinite, since that
761 //    is consistent with the behavior for f == 0, g < 0 and 0 < g < 1.
762 //    Practically any definition could have been justified because mathematical
763 //    consistency has been lost at this point.
764 //
765 // 7. For f < 0, g integer, dg == 0: (f + df)^(g + dg) ~= f^g + g * f^(g - 1) df
766 //    This is equivalent to the case where f is a differentiable function and g
767 //    is a constant (to first order).
768 //
769 // 8. For f < 0, g integer, dg != 0: The value is finite but the derivatives are
770 //    not, because any change in the value of g moves us away from the point
771 //    with a real-valued answer into the region with complex-valued answers.
772 //
773 // 9. For f < 0, g noninteger: The value and derivatives of f^g are not finite.
774
775 template <typename T, int N> inline
776 Jet<T, N> pow(const Jet<T, N>& f, const Jet<T, N>& g) {
777   if (f.a == 0 && g.a >= 1) {
778     // Handle cases 2 and 3.
779     if (g.a > 1) {
780       return Jet<T, N>(T(0.0));
781     }
782     return f;
783   }
784   if (f.a < 0 && g.a == floor(g.a)) {
785     // Handle cases 7 and 8.
786     T const tmp = g.a * pow(f.a, g.a - T(1.0));
787     Jet<T, N> ret(pow(f.a, g.a), tmp * f.v);
788     for (int i = 0; i < N; i++) {
789       if (g.v[i] != T(0.0)) {
790         // Return a NaN when g.v != 0.
791         ret.v[i] = std::numeric_limits<T>::quiet_NaN();
792       }
793     }
794     return ret;
795   }
796   // Handle the remaining cases. For cases 4,5,6,9 we allow the log() function
797   // to generate -HUGE_VAL or NaN, since those cases result in a nonfinite
798   // derivative.
799   T const tmp1 = pow(f.a, g.a);
800   T const tmp2 = g.a * pow(f.a, g.a - T(1.0));
801   T const tmp3 = tmp1 * log(f.a);
802   return Jet<T, N>(tmp1, tmp2 * f.v + tmp3 * g.v);
803 }
804
805 // Define the helper functions Eigen needs to embed Jet types.
806 //
807 // NOTE(keir): machine_epsilon() and precision() are missing, because they don't
808 // work with nested template types (e.g. where the scalar is itself templated).
809 // Among other things, this means that decompositions of Jet's does not work,
810 // for example
811 //
812 //   Matrix<Jet<T, N> ... > A, x, b;
813 //   ...
814 //   A.solve(b, &x)
815 //
816 // does not work and will fail with a strange compiler error.
817 //
818 // TODO(keir): This is an Eigen 2.0 limitation that is lifted in 3.0. When we
819 // switch to 3.0, also add the rest of the specialization functionality.
820 template<typename T, int N> inline const Jet<T, N>& ei_conj(const Jet<T, N>& x) { return x;              }  // NOLINT
821 template<typename T, int N> inline const Jet<T, N>& ei_real(const Jet<T, N>& x) { return x;              }  // NOLINT
822 template<typename T, int N> inline       Jet<T, N>  ei_imag(const Jet<T, N>&  ) { return Jet<T, N>(0.0); }  // NOLINT
823 template<typename T, int N> inline       Jet<T, N>  ei_abs (const Jet<T, N>& x) { return fabs(x);        }  // NOLINT
824 template<typename T, int N> inline       Jet<T, N>  ei_abs2(const Jet<T, N>& x) { return x * x;          }  // NOLINT
825 template<typename T, int N> inline       Jet<T, N>  ei_sqrt(const Jet<T, N>& x) { return sqrt(x);        }  // NOLINT
826 template<typename T, int N> inline       Jet<T, N>  ei_exp (const Jet<T, N>& x) { return exp(x);         }  // NOLINT
827 template<typename T, int N> inline       Jet<T, N>  ei_log (const Jet<T, N>& x) { return log(x);         }  // NOLINT
828 template<typename T, int N> inline       Jet<T, N>  ei_sin (const Jet<T, N>& x) { return sin(x);         }  // NOLINT
829 template<typename T, int N> inline       Jet<T, N>  ei_cos (const Jet<T, N>& x) { return cos(x);         }  // NOLINT
830 template<typename T, int N> inline       Jet<T, N>  ei_tan (const Jet<T, N>& x) { return tan(x);         }  // NOLINT
831 template<typename T, int N> inline       Jet<T, N>  ei_atan(const Jet<T, N>& x) { return atan(x);        }  // NOLINT
832 template<typename T, int N> inline       Jet<T, N>  ei_sinh(const Jet<T, N>& x) { return sinh(x);        }  // NOLINT
833 template<typename T, int N> inline       Jet<T, N>  ei_cosh(const Jet<T, N>& x) { return cosh(x);        }  // NOLINT
834 template<typename T, int N> inline       Jet<T, N>  ei_tanh(const Jet<T, N>& x) { return tanh(x);        }  // NOLINT
835 template<typename T, int N> inline       Jet<T, N>  ei_pow (const Jet<T, N>& x, Jet<T, N> y) { return pow(x, y); }  // NOLINT
836
837 // Note: This has to be in the ceres namespace for argument dependent lookup to
838 // function correctly. Otherwise statements like CHECK_LE(x, 2.0) fail with
839 // strange compile errors.
840 template <typename T, int N>
841 inline std::ostream &operator<<(std::ostream &s, const Jet<T, N>& z) {
842   s << "[" << z.a << " ; ";
843   for (int i = 0; i < N; ++i) {
844     s << z.v[i];
845     if (i != N - 1) {
846       s << ", ";
847     }
848   }
849   s << "]";
850   return s;
851 }
852
853 }  // namespace ceres
854
855 namespace Eigen {
856
857 // Creating a specialization of NumTraits enables placing Jet objects inside
858 // Eigen arrays, getting all the goodness of Eigen combined with autodiff.
859 template<typename T, int N>
860 struct NumTraits<ceres::Jet<T, N> > {
861   typedef ceres::Jet<T, N> Real;
862   typedef ceres::Jet<T, N> NonInteger;
863   typedef ceres::Jet<T, N> Nested;
864   typedef ceres::Jet<T, N> Literal;
865
866   static typename ceres::Jet<T, N> dummy_precision() {
867     return ceres::Jet<T, N>(1e-12);
868   }
869
870   static inline Real epsilon() {
871     return Real(std::numeric_limits<T>::epsilon());
872   }
873
874   static inline int digits10() { return NumTraits<T>::digits10(); }
875
876   enum {
877     IsComplex = 0,
878     IsInteger = 0,
879     IsSigned,
880     ReadCost = 1,
881     AddCost = 1,
882     // For Jet types, multiplication is more expensive than addition.
883     MulCost = 3,
884     HasFloatingPoint = 1,
885     RequireInitialization = 1
886   };
887
888   template<bool Vectorized>
889   struct Div {
890     enum {
891 #if defined(EIGEN_VECTORIZE_AVX)
892       AVX = true,
893 #else
894       AVX = false,
895 #endif
896
897       // Assuming that for Jets, division is as expensive as
898       // multiplication.
899       Cost = 3
900     };
901   };
902 };
903
904 #if EIGEN_VERSION_AT_LEAST(3, 3, 0)
905 // Specifying the return type of binary operations between Jets and scalar types
906 // allows you to perform matrix/array operations with Eigen matrices and arrays
907 // such as addition, subtraction, multiplication, and division where one Eigen
908 // matrix/array is of type Jet and the other is a scalar type. This improves
909 // performance by using the optimized scalar-to-Jet binary operations but
910 // is only available on Eigen versions >= 3.3
911 template <typename BinaryOp, typename T, int N>
912 struct ScalarBinaryOpTraits<ceres::Jet<T, N>, T, BinaryOp> {
913   typedef ceres::Jet<T, N> ReturnType;
914 };
915 template <typename BinaryOp, typename T, int N>
916 struct ScalarBinaryOpTraits<T, ceres::Jet<T, N>, BinaryOp> {
917   typedef ceres::Jet<T, N> ReturnType;
918 };
919 #endif  // EIGEN_VERSION_AT_LEAST(3, 3, 0)
920
921 }  // namespace Eigen
922
923 #endif  // CERES_PUBLIC_JET_H_