blob: 6ae6ef7e29a7f711c8d2eebe557e9e01aa37a75c [file] [log] [blame]
Austin Schuh70cc9552019-01-21 19:46:48 -08001// 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#include "ceres/jet.h"
32
33#include <Eigen/Dense>
34#include <algorithm>
35#include <cmath>
36
37#include "ceres/stringprintf.h"
38#include "ceres/test_util.h"
39#include "glog/logging.h"
40#include "gtest/gtest.h"
41
42#define VL VLOG(1)
43
44namespace ceres {
45namespace internal {
46
47const double kE = 2.71828182845904523536;
48
49typedef Jet<double, 2> J;
50
51// Convenient shorthand for making a jet.
52J MakeJet(double a, double v0, double v1) {
53 J z;
54 z.a = a;
55 z.v[0] = v0;
56 z.v[1] = v1;
57 return z;
58}
59
60// On a 32-bit optimized build, the mismatch is about 1.4e-14.
61double const kTolerance = 1e-13;
62
63void ExpectJetsClose(const J &x, const J &y) {
64 ExpectClose(x.a, y.a, kTolerance);
65 ExpectClose(x.v[0], y.v[0], kTolerance);
66 ExpectClose(x.v[1], y.v[1], kTolerance);
67}
68
69const double kStep = 1e-8;
70const double kNumericalTolerance = 1e-6; // Numeric derivation is quite inexact
71
72// Differentiate using Jet and confirm results with numerical derivation.
73template<typename Function>
74void NumericalTest(const char* name, const Function& f, const double x) {
75 const double exact_dx = f(MakeJet(x, 1.0, 0.0)).v[0];
76 const double estimated_dx =
77 (f(J(x + kStep)).a - f(J(x - kStep)).a) / (2.0 * kStep);
78 VL << name << "(" << x << "), exact dx: "
79 << exact_dx << ", estimated dx: " << estimated_dx;
80 ExpectClose(exact_dx, estimated_dx, kNumericalTolerance);
81}
82
83// Same as NumericalTest, but given a function taking two arguments.
84template<typename Function>
85void NumericalTest2(const char* name, const Function& f,
86 const double x, const double y) {
87 const J exact_delta = f(MakeJet(x, 1.0, 0.0), MakeJet(y, 0.0, 1.0));
88 const double exact_dx = exact_delta.v[0];
89 const double exact_dy = exact_delta.v[1];
90
91 // Sanity check – these should be equivalent:
92 EXPECT_EQ(exact_dx, f(MakeJet(x, 1.0, 0.0), MakeJet(y, 0.0, 0.0)).v[0]);
93 EXPECT_EQ(exact_dx, f(MakeJet(x, 0.0, 1.0), MakeJet(y, 0.0, 0.0)).v[1]);
94 EXPECT_EQ(exact_dy, f(MakeJet(x, 0.0, 0.0), MakeJet(y, 1.0, 0.0)).v[0]);
95 EXPECT_EQ(exact_dy, f(MakeJet(x, 0.0, 0.0), MakeJet(y, 0.0, 1.0)).v[1]);
96
97 const double estimated_dx =
98 (f(J(x + kStep), J(y)).a - f(J(x - kStep), J(y)).a) / (2.0 * kStep);
99 const double estimated_dy =
100 (f(J(x), J(y + kStep)).a - f(J(x), J(y - kStep)).a) / (2.0 * kStep);
101 VL << name << "(" << x << ", " << y << "), exact dx: "
102 << exact_dx << ", estimated dx: " << estimated_dx;
103 ExpectClose(exact_dx, estimated_dx, kNumericalTolerance);
104 VL << name << "(" << x << ", " << y << "), exact dy: "
105 << exact_dy << ", estimated dy: " << estimated_dy;
106 ExpectClose(exact_dy, estimated_dy, kNumericalTolerance);
107}
108
109TEST(Jet, Jet) {
110 // Pick arbitrary values for x and y.
111 J x = MakeJet(2.3, -2.7, 1e-3);
112 J y = MakeJet(1.7, 0.5, 1e+2);
113
114 VL << "x = " << x;
115 VL << "y = " << y;
116
117 { // Check that log(exp(x)) == x.
118 J z = exp(x);
119 J w = log(z);
120 VL << "z = " << z;
121 VL << "w = " << w;
122 ExpectJetsClose(w, x);
123 }
124
125 { // Check that (x * y) / x == y.
126 J z = x * y;
127 J w = z / x;
128 VL << "z = " << z;
129 VL << "w = " << w;
130 ExpectJetsClose(w, y);
131 }
132
133 { // Check that sqrt(x * x) == x.
134 J z = x * x;
135 J w = sqrt(z);
136 VL << "z = " << z;
137 VL << "w = " << w;
138 ExpectJetsClose(w, x);
139 }
140
141 { // Check that sqrt(y) * sqrt(y) == y.
142 J z = sqrt(y);
143 J w = z * z;
144 VL << "z = " << z;
145 VL << "w = " << w;
146 ExpectJetsClose(w, y);
147 }
148
149 NumericalTest("sqrt", sqrt<double, 2>, 0.00001);
150 NumericalTest("sqrt", sqrt<double, 2>, 1.0);
151
152 { // Check that cos(2*x) = cos(x)^2 - sin(x)^2
153 J z = cos(J(2.0) * x);
154 J w = cos(x)*cos(x) - sin(x)*sin(x);
155 VL << "z = " << z;
156 VL << "w = " << w;
157 ExpectJetsClose(w, z);
158 }
159
160 { // Check that sin(2*x) = 2*cos(x)*sin(x)
161 J z = sin(J(2.0) * x);
162 J w = J(2.0)*cos(x)*sin(x);
163 VL << "z = " << z;
164 VL << "w = " << w;
165 ExpectJetsClose(w, z);
166 }
167
168 { // Check that cos(x)*cos(x) + sin(x)*sin(x) = 1
169 J z = cos(x) * cos(x);
170 J w = sin(x) * sin(x);
171 VL << "z = " << z;
172 VL << "w = " << w;
173 ExpectJetsClose(z + w, J(1.0));
174 }
175
176 { // Check that atan2(r*sin(t), r*cos(t)) = t.
177 J t = MakeJet(0.7, -0.3, +1.5);
178 J r = MakeJet(2.3, 0.13, -2.4);
179 VL << "t = " << t;
180 VL << "r = " << r;
181
182 J u = atan2(r * sin(t), r * cos(t));
183 VL << "u = " << u;
184
185 ExpectJetsClose(u, t);
186 }
187
188 { // Check that tan(x) = sin(x) / cos(x).
189 J z = tan(x);
190 J w = sin(x) / cos(x);
191 VL << "z = " << z;
192 VL << "w = " << w;
193 ExpectJetsClose(z, w);
194 }
195
196 { // Check that tan(atan(x)) = x.
197 J z = tan(atan(x));
198 J w = x;
199 VL << "z = " << z;
200 VL << "w = " << w;
201 ExpectJetsClose(z, w);
202 }
203
204 { // Check that cosh(x)*cosh(x) - sinh(x)*sinh(x) = 1
205 J z = cosh(x) * cosh(x);
206 J w = sinh(x) * sinh(x);
207 VL << "z = " << z;
208 VL << "w = " << w;
209 ExpectJetsClose(z - w, J(1.0));
210 }
211
212 { // Check that tanh(x + y) = (tanh(x) + tanh(y)) / (1 + tanh(x) tanh(y))
213 J z = tanh(x + y);
214 J w = (tanh(x) + tanh(y)) / (J(1.0) + tanh(x) * tanh(y));
215 VL << "z = " << z;
216 VL << "w = " << w;
217 ExpectJetsClose(z, w);
218 }
219
220 { // Check that pow(x, 1) == x.
221 VL << "x = " << x;
222
223 J u = pow(x, 1.);
224 VL << "u = " << u;
225
226 ExpectJetsClose(x, u);
227 }
228
229 { // Check that pow(x, 1) == x.
230 J y = MakeJet(1, 0.0, 0.0);
231 VL << "x = " << x;
232 VL << "y = " << y;
233
234 J u = pow(x, y);
235 VL << "u = " << u;
236
237 ExpectJetsClose(x, u);
238 }
239
240 { // Check that pow(e, log(x)) == x.
241 J logx = log(x);
242
243 VL << "x = " << x;
244 VL << "y = " << y;
245
246 J u = pow(kE, logx);
247 VL << "u = " << u;
248
249 ExpectJetsClose(x, u);
250 }
251
252 { // Check that pow(e, log(x)) == x.
253 J logx = log(x);
254 J e = MakeJet(kE, 0., 0.);
255 VL << "x = " << x;
256 VL << "log(x) = " << logx;
257
258 J u = pow(e, logx);
259 VL << "u = " << u;
260
261 ExpectJetsClose(x, u);
262 }
263
264 { // Check that pow(e, log(x)) == x.
265 J logx = log(x);
266 J e = MakeJet(kE, 0., 0.);
267 VL << "x = " << x;
268 VL << "logx = " << logx;
269
270 J u = pow(e, logx);
271 VL << "u = " << u;
272
273 ExpectJetsClose(x, u);
274 }
275
276 { // Check that pow(x,y) = exp(y*log(x)).
277 J logx = log(x);
278 J e = MakeJet(kE, 0., 0.);
279 VL << "x = " << x;
280 VL << "logx = " << logx;
281
282 J u = pow(e, y*logx);
283 J v = pow(x, y);
284 VL << "u = " << u;
285 VL << "v = " << v;
286
287 ExpectJetsClose(v, u);
288 }
289
290 { // Check that pow(0, y) == 0 for y > 1, with both arguments Jets.
291 // This tests special case handling inside pow().
292 J a = MakeJet(0, 1, 2);
293 J b = MakeJet(2, 3, 4);
294 VL << "a = " << a;
295 VL << "b = " << b;
296
297 J c = pow(a, b);
298 VL << "a^b = " << c;
299 ExpectJetsClose(c, MakeJet(0, 0, 0));
300 }
301
302 { // Check that pow(0, y) == 0 for y == 1, with both arguments Jets.
303 // This tests special case handling inside pow().
304 J a = MakeJet(0, 1, 2);
305 J b = MakeJet(1, 3, 4);
306 VL << "a = " << a;
307 VL << "b = " << b;
308
309 J c = pow(a, b);
310 VL << "a^b = " << c;
311 ExpectJetsClose(c, MakeJet(0, 1, 2));
312 }
313
314 { // Check that pow(0, <1) is not finite, with both arguments Jets.
315 for (int i = 1; i < 10; i++) {
316 J a = MakeJet(0, 1, 2);
317 J b = MakeJet(i*0.1, 3, 4); // b = 0.1 ... 0.9
318 VL << "a = " << a;
319 VL << "b = " << b;
320
321 J c = pow(a, b);
322 VL << "a^b = " << c;
323 EXPECT_EQ(c.a, 0.0);
324 EXPECT_FALSE(IsFinite(c.v[0]));
325 EXPECT_FALSE(IsFinite(c.v[1]));
326 }
327 for (int i = -10; i < 0; i++) {
328 J a = MakeJet(0, 1, 2);
329 J b = MakeJet(i*0.1, 3, 4); // b = -1,-0.9 ... -0.1
330 VL << "a = " << a;
331 VL << "b = " << b;
332
333 J c = pow(a, b);
334 VL << "a^b = " << c;
335 EXPECT_FALSE(IsFinite(c.a));
336 EXPECT_FALSE(IsFinite(c.v[0]));
337 EXPECT_FALSE(IsFinite(c.v[1]));
338 }
339
340 {
341 // The special case of 0^0 = 1 defined by the C standard.
342 J a = MakeJet(0, 1, 2);
343 J b = MakeJet(0, 3, 4);
344 VL << "a = " << a;
345 VL << "b = " << b;
346
347 J c = pow(a, b);
348 VL << "a^b = " << c;
349 EXPECT_EQ(c.a, 1.0);
350 EXPECT_FALSE(IsFinite(c.v[0]));
351 EXPECT_FALSE(IsFinite(c.v[1]));
352 }
353 }
354
355 { // Check that pow(<0, b) is correct for integer b.
356 // This tests special case handling inside pow().
357 J a = MakeJet(-1.5, 3, 4);
358
359 // b integer:
360 for (int i = -10; i <= 10; i++) {
361 J b = MakeJet(i, 0, 5);
362 VL << "a = " << a;
363 VL << "b = " << b;
364
365 J c = pow(a, b);
366 VL << "a^b = " << c;
367 ExpectClose(c.a, pow(-1.5, i), kTolerance);
368 EXPECT_TRUE(IsFinite(c.v[0]));
369 EXPECT_FALSE(IsFinite(c.v[1]));
370 ExpectClose(c.v[0], i * pow(-1.5, i - 1) * 3.0, kTolerance);
371 }
372 }
373
374 { // Check that pow(<0, b) is correct for noninteger b.
375 // This tests special case handling inside pow().
376 J a = MakeJet(-1.5, 3, 4);
377 J b = MakeJet(-2.5, 0, 5);
378 VL << "a = " << a;
379 VL << "b = " << b;
380
381 J c = pow(a, b);
382 VL << "a^b = " << c;
383 EXPECT_FALSE(IsFinite(c.a));
384 EXPECT_FALSE(IsFinite(c.v[0]));
385 EXPECT_FALSE(IsFinite(c.v[1]));
386 }
387
388 {
389 // Check that pow(0,y) == 0 for y == 2, with the second argument a
390 // Jet. This tests special case handling inside pow().
391 double a = 0;
392 J b = MakeJet(2, 3, 4);
393 VL << "a = " << a;
394 VL << "b = " << b;
395
396 J c = pow(a, b);
397 VL << "a^b = " << c;
398 ExpectJetsClose(c, MakeJet(0, 0, 0));
399 }
400
401 {
402 // Check that pow(<0,y) is correct for integer y. This tests special case
403 // handling inside pow().
404 double a = -1.5;
405 for (int i = -10; i <= 10; i++) {
406 J b = MakeJet(i, 3, 0);
407 VL << "a = " << a;
408 VL << "b = " << b;
409
410 J c = pow(a, b);
411 VL << "a^b = " << c;
412 ExpectClose(c.a, pow(-1.5, i), kTolerance);
413 EXPECT_FALSE(IsFinite(c.v[0]));
414 EXPECT_TRUE(IsFinite(c.v[1]));
415 ExpectClose(c.v[1], 0, kTolerance);
416 }
417 }
418
419 {
420 // Check that pow(<0,y) is correct for noninteger y. This tests special
421 // case handling inside pow().
422 double a = -1.5;
423 J b = MakeJet(-3.14, 3, 0);
424 VL << "a = " << a;
425 VL << "b = " << b;
426
427 J c = pow(a, b);
428 VL << "a^b = " << c;
429 EXPECT_FALSE(IsFinite(c.a));
430 EXPECT_FALSE(IsFinite(c.v[0]));
431 EXPECT_FALSE(IsFinite(c.v[1]));
432 }
433
434 { // Check that 1 + x == x + 1.
435 J a = x + 1.0;
436 J b = 1.0 + x;
437 J c = x;
438 c += 1.0;
439
440 ExpectJetsClose(a, b);
441 ExpectJetsClose(a, c);
442 }
443
444 { // Check that 1 - x == -(x - 1).
445 J a = 1.0 - x;
446 J b = -(x - 1.0);
447 J c = x;
448 c -= 1.0;
449
450 ExpectJetsClose(a, b);
451 ExpectJetsClose(a, -c);
452 }
453
454 { // Check that (x/s)*s == (x*s)/s.
455 J a = x / 5.0;
456 J b = x * 5.0;
457 J c = x;
458 c /= 5.0;
459 J d = x;
460 d *= 5.0;
461
462 ExpectJetsClose(5.0 * a, b / 5.0);
463 ExpectJetsClose(a, c);
464 ExpectJetsClose(b, d);
465 }
466
467 { // Check that x / y == 1 / (y / x).
468 J a = x / y;
469 J b = 1.0 / (y / x);
470 VL << "a = " << a;
471 VL << "b = " << b;
472
473 ExpectJetsClose(a, b);
474 }
475
476 { // Check that abs(-x * x) == sqrt(x * x).
477 ExpectJetsClose(abs(-x), sqrt(x * x));
478 }
479
480 { // Check that cos(acos(x)) == x.
481 J a = MakeJet(0.1, -2.7, 1e-3);
482 ExpectJetsClose(cos(acos(a)), a);
483 ExpectJetsClose(acos(cos(a)), a);
484
485 J b = MakeJet(0.6, 0.5, 1e+2);
486 ExpectJetsClose(cos(acos(b)), b);
487 ExpectJetsClose(acos(cos(b)), b);
488 }
489
490 { // Check that sin(asin(x)) == x.
491 J a = MakeJet(0.1, -2.7, 1e-3);
492 ExpectJetsClose(sin(asin(a)), a);
493 ExpectJetsClose(asin(sin(a)), a);
494
495 J b = MakeJet(0.4, 0.5, 1e+2);
496 ExpectJetsClose(sin(asin(b)), b);
497 ExpectJetsClose(asin(sin(b)), b);
498 }
499
500 {
501 J zero = J(0.0);
502
503 // Check that J0(0) == 1.
504 ExpectJetsClose(BesselJ0(zero), J(1.0));
505
506 // Check that J1(0) == 0.
507 ExpectJetsClose(BesselJ1(zero), zero);
508
509 // Check that J2(0) == 0.
510 ExpectJetsClose(BesselJn(2, zero), zero);
511
512 // Check that J3(0) == 0.
513 ExpectJetsClose(BesselJn(3, zero), zero);
514
515 J z = MakeJet(0.1, -2.7, 1e-3);
516
517 // Check that J0(z) == Jn(0,z).
518 ExpectJetsClose(BesselJ0(z), BesselJn(0, z));
519
520 // Check that J1(z) == Jn(1,z).
521 ExpectJetsClose(BesselJ1(z), BesselJn(1, z));
522
523 // Check that J0(z)+J2(z) == (2/z)*J1(z).
524 // See formula http://dlmf.nist.gov/10.6.E1
525 ExpectJetsClose(BesselJ0(z) + BesselJn(2, z), (2.0 / z) * BesselJ1(z));
526 }
527
528 { // Check that floor of a positive number works.
529 J a = MakeJet(0.1, -2.7, 1e-3);
530 J b = floor(a);
531 J expected = MakeJet(floor(a.a), 0.0, 0.0);
532 EXPECT_EQ(expected, b);
533 }
534
535 { // Check that floor of a negative number works.
536 J a = MakeJet(-1.1, -2.7, 1e-3);
537 J b = floor(a);
538 J expected = MakeJet(floor(a.a), 0.0, 0.0);
539 EXPECT_EQ(expected, b);
540 }
541
542 { // Check that floor of a positive number works.
543 J a = MakeJet(10.123, -2.7, 1e-3);
544 J b = floor(a);
545 J expected = MakeJet(floor(a.a), 0.0, 0.0);
546 EXPECT_EQ(expected, b);
547 }
548
549 { // Check that ceil of a positive number works.
550 J a = MakeJet(0.1, -2.7, 1e-3);
551 J b = ceil(a);
552 J expected = MakeJet(ceil(a.a), 0.0, 0.0);
553 EXPECT_EQ(expected, b);
554 }
555
556 { // Check that ceil of a negative number works.
557 J a = MakeJet(-1.1, -2.7, 1e-3);
558 J b = ceil(a);
559 J expected = MakeJet(ceil(a.a), 0.0, 0.0);
560 EXPECT_EQ(expected, b);
561 }
562
563 { // Check that ceil of a positive number works.
564 J a = MakeJet(10.123, -2.7, 1e-3);
565 J b = ceil(a);
566 J expected = MakeJet(ceil(a.a), 0.0, 0.0);
567 EXPECT_EQ(expected, b);
568 }
569
570 { // Check that cbrt(x * x * x) == x.
571 J z = x * x * x;
572 J w = cbrt(z);
573 VL << "z = " << z;
574 VL << "w = " << w;
575 ExpectJetsClose(w, x);
576 }
577
578 { // Check that cbrt(y) * cbrt(y) * cbrt(y) == y.
579 J z = cbrt(y);
580 J w = z * z * z;
581 VL << "z = " << z;
582 VL << "w = " << w;
583 ExpectJetsClose(w, y);
584 }
585
586 { // Check that cbrt(x) == pow(x, 1/3).
587 J z = cbrt(x);
588 J w = pow(x, 1.0 / 3.0);
589 VL << "z = " << z;
590 VL << "w = " << w;
591 ExpectJetsClose(z, w);
592 }
593 NumericalTest("cbrt", cbrt<double, 2>, -1.0);
594 NumericalTest("cbrt", cbrt<double, 2>, -1e-5);
595 NumericalTest("cbrt", cbrt<double, 2>, 1e-5);
596 NumericalTest("cbrt", cbrt<double, 2>, 1.0);
597
598 { // Check that exp2(x) == exp(x * log(2))
599 J z = exp2(x);
600 J w = exp(x * log(2.0));
601 VL << "z = " << z;
602 VL << "w = " << w;
603 ExpectJetsClose(z, w);
604 }
605 NumericalTest("exp2", exp2<double, 2>, -1.0);
606 NumericalTest("exp2", exp2<double, 2>, -1e-5);
607 NumericalTest("exp2", exp2<double, 2>, -1e-200);
608 NumericalTest("exp2", exp2<double, 2>, 0.0);
609 NumericalTest("exp2", exp2<double, 2>, 1e-200);
610 NumericalTest("exp2", exp2<double, 2>, 1e-5);
611 NumericalTest("exp2", exp2<double, 2>, 1.0);
612
613 { // Check that log2(x) == log(x) / log(2)
614 J z = log2(x);
615 J w = log(x) / log(2.0);
616 VL << "z = " << z;
617 VL << "w = " << w;
618 ExpectJetsClose(z, w);
619 }
620 NumericalTest("log2", log2<double, 2>, 1e-5);
621 NumericalTest("log2", log2<double, 2>, 1.0);
622 NumericalTest("log2", log2<double, 2>, 100.0);
623
624 { // Check that hypot(x, y) == sqrt(x^2 + y^2)
625 J h = hypot(x, y);
626 J s = sqrt(x*x + y*y);
627 VL << "h = " << h;
628 VL << "s = " << s;
629 ExpectJetsClose(h, s);
630 }
631
632 { // Check that hypot(x, x) == sqrt(2) * abs(x)
633 J h = hypot(x, x);
634 J s = sqrt(2.0) * abs(x);
635 VL << "h = " << h;
636 VL << "s = " << s;
637 ExpectJetsClose(h, s);
638 }
639
640 { // Check that the derivative is zero tangentially to the circle:
641 J h = hypot(MakeJet(2.0, 1.0, 1.0), MakeJet(2.0, 1.0, -1.0));
642 VL << "h = " << h;
643 ExpectJetsClose(h, MakeJet(sqrt(8.0), std::sqrt(2.0), 0.0));
644 }
645
646 { // Check that hypot(x, 0) == x
647 J zero = MakeJet(0.0, 2.0, 3.14);
648 J h = hypot(x, zero);
649 VL << "h = " << h;
650 ExpectJetsClose(x, h);
651 }
652
653 { // Check that hypot(0, y) == y
654 J zero = MakeJet(0.0, 2.0, 3.14);
655 J h = hypot(zero, y);
656 VL << "h = " << h;
657 ExpectJetsClose(y, h);
658 }
659
660 { // Check that hypot(x, 0) == sqrt(x * x) == x, even when x * x underflows:
661 EXPECT_EQ(DBL_MIN * DBL_MIN, 0.0); // Make sure it underflows
662 J huge = MakeJet(DBL_MIN, 2.0, 3.14);
663 J h = hypot(huge, J(0.0));
664 VL << "h = " << h;
665 ExpectJetsClose(h, huge);
666 }
667
668 { // Check that hypot(x, 0) == sqrt(x * x) == x, even when x * x overflows:
669 EXPECT_EQ(DBL_MAX * DBL_MAX, std::numeric_limits<double>::infinity());
670 J huge = MakeJet(DBL_MAX, 2.0, 3.14);
671 J h = hypot(huge, J(0.0));
672 VL << "h = " << h;
673 ExpectJetsClose(h, huge);
674 }
675
676 NumericalTest2("hypot", hypot<double, 2>, 0.0, 1e-5);
677 NumericalTest2("hypot", hypot<double, 2>, -1e-5, 0.0);
678 NumericalTest2("hypot", hypot<double, 2>, 1e-5, 1e-5);
679 NumericalTest2("hypot", hypot<double, 2>, 0.0, 1.0);
680 NumericalTest2("hypot", hypot<double, 2>, 1e-3, 1.0);
681 NumericalTest2("hypot", hypot<double, 2>, 1e-3, -1.0);
682 NumericalTest2("hypot", hypot<double, 2>, -1e-3, 1.0);
683 NumericalTest2("hypot", hypot<double, 2>, -1e-3, -1.0);
684 NumericalTest2("hypot", hypot<double, 2>, 1.0, 2.0);
685
686 {
687 J z = fmax(x, y);
688 VL << "z = " << z;
689 ExpectJetsClose(x, z);
690 }
691
692 {
693 J z = fmin(x, y);
694 VL << "z = " << z;
695 ExpectJetsClose(y, z);
696 }
697
698}
699
700TEST(Jet, JetsInEigenMatrices) {
701 J x = MakeJet(2.3, -2.7, 1e-3);
702 J y = MakeJet(1.7, 0.5, 1e+2);
703 J z = MakeJet(5.3, -4.7, 1e-3);
704 J w = MakeJet(9.7, 1.5, 10.1);
705
706 Eigen::Matrix<J, 2, 2> M;
707 Eigen::Matrix<J, 2, 1> v, r1, r2;
708
709 M << x, y, z, w;
710 v << x, z;
711
712 // Check that M * v == (v^T * M^T)^T
713 r1 = M * v;
714 r2 = (v.transpose() * M.transpose()).transpose();
715
716 ExpectJetsClose(r1(0), r2(0));
717 ExpectJetsClose(r1(1), r2(1));
718}
719
720TEST(JetTraitsTest, ClassificationMixed) {
721 Jet<double, 3> a(5.5, 0);
722 a.v[0] = std::numeric_limits<double>::quiet_NaN();
723 a.v[1] = std::numeric_limits<double>::infinity();
724 a.v[2] = -std::numeric_limits<double>::infinity();
725 EXPECT_FALSE(IsFinite(a));
726 EXPECT_FALSE(IsNormal(a));
727 EXPECT_TRUE(IsInfinite(a));
728 EXPECT_TRUE(IsNaN(a));
729}
730
731TEST(JetTraitsTest, ClassificationNaN) {
732 Jet<double, 3> a(5.5, 0);
733 a.v[0] = std::numeric_limits<double>::quiet_NaN();
734 a.v[1] = 0.0;
735 a.v[2] = 0.0;
736 EXPECT_FALSE(IsFinite(a));
737 EXPECT_FALSE(IsNormal(a));
738 EXPECT_FALSE(IsInfinite(a));
739 EXPECT_TRUE(IsNaN(a));
740}
741
742TEST(JetTraitsTest, ClassificationInf) {
743 Jet<double, 3> a(5.5, 0);
744 a.v[0] = std::numeric_limits<double>::infinity();
745 a.v[1] = 0.0;
746 a.v[2] = 0.0;
747 EXPECT_FALSE(IsFinite(a));
748 EXPECT_FALSE(IsNormal(a));
749 EXPECT_TRUE(IsInfinite(a));
750 EXPECT_FALSE(IsNaN(a));
751}
752
753TEST(JetTraitsTest, ClassificationFinite) {
754 Jet<double, 3> a(5.5, 0);
755 a.v[0] = 100.0;
756 a.v[1] = 1.0;
757 a.v[2] = 3.14159;
758 EXPECT_TRUE(IsFinite(a));
759 EXPECT_TRUE(IsNormal(a));
760 EXPECT_FALSE(IsInfinite(a));
761 EXPECT_FALSE(IsNaN(a));
762}
763
764#if EIGEN_VERSION_AT_LEAST(3, 3, 0)
765
766// The following test ensures that Jets have all the appropriate Eigen
767// related traits so that they can be used as part of matrix
768// decompositions.
769TEST(Jet, FullRankEigenLLTSolve) {
770 Eigen::Matrix<J, 3, 3> A;
771 Eigen::Matrix<J, 3, 1> b, x;
772 for (int i = 0; i < 3; ++i) {
773 for (int j = 0; j < 3; ++j) {
774 A(i,j) = MakeJet(0.0, i, j * j);
775 }
776 b(i) = MakeJet(i, i, i);
777 x(i) = MakeJet(0.0, 0.0, 0.0);
778 A(i,i) = MakeJet(1.0, i, i * i);
779 }
780 x = A.llt().solve(b);
781 for (int i = 0; i < 3; ++i) {
782 EXPECT_EQ(x(i).a, b(i).a);
783 }
784}
785
786TEST(Jet, FullRankEigenLDLTSolve) {
787 Eigen::Matrix<J, 3, 3> A;
788 Eigen::Matrix<J, 3, 1> b, x;
789 for (int i = 0; i < 3; ++i) {
790 for (int j = 0; j < 3; ++j) {
791 A(i,j) = MakeJet(0.0, i, j * j);
792 }
793 b(i) = MakeJet(i, i, i);
794 x(i) = MakeJet(0.0, 0.0, 0.0);
795 A(i,i) = MakeJet(1.0, i, i * i);
796 }
797 x = A.ldlt().solve(b);
798 for (int i = 0; i < 3; ++i) {
799 EXPECT_EQ(x(i).a, b(i).a);
800 }
801}
802
803TEST(Jet, FullRankEigenLUSolve) {
804 Eigen::Matrix<J, 3, 3> A;
805 Eigen::Matrix<J, 3, 1> b, x;
806 for (int i = 0; i < 3; ++i) {
807 for (int j = 0; j < 3; ++j) {
808 A(i,j) = MakeJet(0.0, i, j * j);
809 }
810 b(i) = MakeJet(i, i, i);
811 x(i) = MakeJet(0.0, 0.0, 0.0);
812 A(i,i) = MakeJet(1.0, i, i * i);
813 }
814
815 x = A.lu().solve(b);
816 for (int i = 0; i < 3; ++i) {
817 EXPECT_EQ(x(i).a, b(i).a);
818 }
819}
820
821// ScalarBinaryOpTraits is only supported on Eigen versions >= 3.3
822TEST(JetTraitsTest, MatrixScalarUnaryOps) {
823 const J x = MakeJet(2.3, -2.7, 1e-3);
824 const J y = MakeJet(1.7, 0.5, 1e+2);
825 Eigen::Matrix<J, 2, 1> a;
826 a << x, y;
827
828 const J sum = a.sum();
829 const J sum2 = a(0) + a(1);
830 ExpectJetsClose(sum, sum2);
831}
832
833TEST(JetTraitsTest, MatrixScalarBinaryOps) {
834 const J x = MakeJet(2.3, -2.7, 1e-3);
835 const J y = MakeJet(1.7, 0.5, 1e+2);
836 const J z = MakeJet(5.3, -4.7, 1e-3);
837 const J w = MakeJet(9.7, 1.5, 10.1);
838
839 Eigen::Matrix<J, 2, 2> M;
840 Eigen::Vector2d v;
841
842 M << x, y, z, w;
843 v << 0.6, -2.1;
844
845 // Check that M * v == M * v.cast<J>().
846 const Eigen::Matrix<J, 2, 1> r1 = M * v;
847 const Eigen::Matrix<J, 2, 1> r2 = M * v.cast<J>();
848
849 ExpectJetsClose(r1(0), r2(0));
850 ExpectJetsClose(r1(1), r2(1));
851
852 // Check that M * a == M * T(a).
853 const double a = 3.1;
854 const Eigen::Matrix<J, 2, 2> r3 = M * a;
855 const Eigen::Matrix<J, 2, 2> r4 = M * J(a);
856
857 ExpectJetsClose(r3(0, 0), r4(0, 0));
858 ExpectJetsClose(r3(1, 0), r4(1, 0));
859 ExpectJetsClose(r3(0, 1), r4(0, 1));
860 ExpectJetsClose(r3(1, 1), r4(1, 1));
861}
862
863TEST(JetTraitsTest, ArrayScalarUnaryOps) {
864 const J x = MakeJet(2.3, -2.7, 1e-3);
865 const J y = MakeJet(1.7, 0.5, 1e+2);
866 Eigen::Array<J, 2, 1> a;
867 a << x, y;
868
869 const J sum = a.sum();
870 const J sum2 = a(0) + a(1);
871 ExpectJetsClose(sum, sum2);
872}
873
874TEST(JetTraitsTest, ArrayScalarBinaryOps) {
875 const J x = MakeJet(2.3, -2.7, 1e-3);
876 const J y = MakeJet(1.7, 0.5, 1e+2);
877
878 Eigen::Array<J, 2, 1> a;
879 Eigen::Array2d b;
880
881 a << x, y;
882 b << 0.6, -2.1;
883
884 // Check that a * b == a * b.cast<T>()
885 const Eigen::Array<J, 2, 1> r1 = a * b;
886 const Eigen::Array<J, 2, 1> r2 = a * b.cast<J>();
887
888 ExpectJetsClose(r1(0), r2(0));
889 ExpectJetsClose(r1(1), r2(1));
890
891 // Check that a * c == a * T(c).
892 const double c = 3.1;
893 const Eigen::Array<J, 2, 1> r3 = a * c;
894 const Eigen::Array<J, 2, 1> r4 = a * J(c);
895
896 ExpectJetsClose(r3(0), r3(0));
897 ExpectJetsClose(r4(1), r4(1));
898}
899#endif // EIGEN_VERSION_AT_LEAST(3, 3, 0)
900
901} // namespace internal
902} // namespace ceres