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