Simbody
3.4 (development)
|
00001 #ifndef SimTK_SimTKCOMMON_TESTING_H_ 00002 #define SimTK_SimTKCOMMON_TESTING_H_ 00003 00004 /* -------------------------------------------------------------------------- * 00005 * Simbody(tm): SimTKcommon * 00006 * -------------------------------------------------------------------------- * 00007 * This is part of the SimTK biosimulation toolkit originating from * 00008 * Simbios, the NIH National Center for Physics-Based Simulation of * 00009 * Biological Structures at Stanford, funded under the NIH Roadmap for * 00010 * Medical Research, grant U54 GM072970. See https://simtk.org/home/simbody. * 00011 * * 00012 * Portions copyright (c) 2009-12 Stanford University and the Authors. * 00013 * Authors: Michael Sherman * 00014 * Contributors: * 00015 * * 00016 * Licensed under the Apache License, Version 2.0 (the "License"); you may * 00017 * not use this file except in compliance with the License. You may obtain a * 00018 * copy of the License at http://www.apache.org/licenses/LICENSE-2.0. * 00019 * * 00020 * Unless required by applicable law or agreed to in writing, software * 00021 * distributed under the License is distributed on an "AS IS" BASIS, * 00022 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied. * 00023 * See the License for the specific language governing permissions and * 00024 * limitations under the License. * 00025 * -------------------------------------------------------------------------- */ 00026 00027 #include "SimTKcommon/basics.h" 00028 #include "SimTKcommon/Simmatrix.h" 00029 #include "SimTKcommon/internal/Random.h" 00030 00031 #include <cmath> 00032 #include <ctime> 00033 #include <algorithm> 00034 #include <iostream> 00035 00041 namespace SimTK { 00042 00155 00156 00157 00158 00159 class Test { 00160 public: 00161 class Subtest; 00162 Test(const std::string& name) : testName(name) 00163 { 00164 std::clog << "Starting test " << testName << " ...\n"; 00165 startTime = std::clock(); 00166 } 00167 ~Test() { 00168 std::clog << "Done. " << testName << " time: " 00169 << 1000*(std::clock()-startTime)/CLOCKS_PER_SEC << "ms.\n"; 00170 } 00171 00172 template <class T> 00173 static double defTol() {return (double)NTraits<typename CNT<T>::Precision>::getSignificant();} 00174 00175 // For dissimilar types, the default tolerance is the narrowest of the two. 00176 template <class T1, class T2> 00177 static double defTol2() {return std::max(defTol<T1>(), defTol<T2>());} 00178 00179 // Scale by the magnitude of the quantities being compared, so that we don't 00180 // ask for unreasonable precision. For magnitudes near zero, we'll be satisfied 00181 // if both are very small without demanding that they must also be relatively 00182 // close. That is, we use a relative tolerance for big numbers and an absolute 00183 // tolerance for small ones. 00184 static bool numericallyEqual(float v1, float v2, int n, double tol=defTol<float>()) { 00185 const float scale = n*std::max(std::max(std::abs(v1), std::abs(v2)), 1.0f); 00186 return std::abs(v1-v2) < scale*(float)tol; 00187 } 00188 static bool numericallyEqual(double v1, double v2, int n, double tol=defTol<double>()) { 00189 const double scale = n*std::max(std::max(std::abs(v1), std::abs(v2)), 1.0); 00190 return std::abs(v1-v2) < scale*(double)tol; 00191 } 00192 static bool numericallyEqual(long double v1, long double v2, int n, double tol=defTol<long double>()) { 00193 const long double scale = n*std::max(std::max(std::abs(v1), std::abs(v2)), 1.0l); 00194 return std::abs(v1-v2) < scale*(long double)tol; 00195 } 00196 00197 // For integers we ignore tolerance. 00198 static bool numericallyEqual(int i1, int i2, int n, double tol=0) {return i1==i2;} 00199 static bool numericallyEqual(unsigned u1, unsigned u2, int n, double tol=0) {return u1==u2;} 00200 00201 // Mixed floating types use default tolerance for the narrower type. 00202 static bool numericallyEqual(float v1, double v2, int n, double tol=defTol<float>()) 00203 { return numericallyEqual((double)v1, v2, n, tol); } 00204 static bool numericallyEqual(double v1, float v2, int n, double tol=defTol<float>()) 00205 { return numericallyEqual(v1, (double)v2, n, tol); } 00206 static bool numericallyEqual(float v1, long double v2, int n, double tol=defTol<float>()) 00207 { return numericallyEqual((long double)v1, v2, n, tol); } 00208 static bool numericallyEqual(long double v1, float v2, int n, double tol=defTol<float>()) 00209 { return numericallyEqual(v1, (long double)v2, n, tol); } 00210 static bool numericallyEqual(double v1, long double v2, int n, double tol=defTol<double>()) 00211 { return numericallyEqual((long double)v1, v2, n, tol); } 00212 static bool numericallyEqual(long double v1, double v2, int n, double tol=defTol<double>()) 00213 { return numericallyEqual(v1, (long double)v2, n, tol); } 00214 00215 // Mixed int/floating just upgrades int to floating type. 00216 static bool numericallyEqual(int i1, float f2, int n, double tol=defTol<float>()) 00217 { return numericallyEqual((float)i1,f2,n,tol); } 00218 static bool numericallyEqual(float f1, int i2, int n, double tol=defTol<float>()) 00219 { return numericallyEqual(f1,(float)i2,n,tol); } 00220 static bool numericallyEqual(unsigned i1, float f2, int n, double tol=defTol<float>()) 00221 { return numericallyEqual((float)i1,f2,n,tol); } 00222 static bool numericallyEqual(float f1, unsigned i2, int n, double tol=defTol<float>()) 00223 { return numericallyEqual(f1,(float)i2,n,tol); } 00224 static bool numericallyEqual(int i1, double f2, int n, double tol=defTol<double>()) 00225 { return numericallyEqual((double)i1,f2,n,tol); } 00226 static bool numericallyEqual(double f1, int i2, int n, double tol=defTol<double>()) 00227 { return numericallyEqual(f1,(double)i2,n,tol); } 00228 static bool numericallyEqual(unsigned i1, double f2, int n, double tol=defTol<double>()) 00229 { return numericallyEqual((double)i1,f2,n,tol); } 00230 static bool numericallyEqual(double f1, unsigned i2, int n, double tol=defTol<double>()) 00231 { return numericallyEqual(f1,(double)i2,n,tol); } 00232 static bool numericallyEqual(int i1, long double f2, int n, double tol=defTol<long double>()) 00233 { return numericallyEqual((long double)i1,f2,n,tol); } 00234 static bool numericallyEqual(long double f1, int i2, int n, double tol=defTol<long double>()) 00235 { return numericallyEqual(f1,(long double)i2,n,tol); } 00236 static bool numericallyEqual(unsigned i1, long double f2, int n, double tol=defTol<long double>()) 00237 { return numericallyEqual((long double)i1,f2,n,tol); } 00238 static bool numericallyEqual(long double f1, unsigned i2, int n, double tol=defTol<long double>()) 00239 { return numericallyEqual(f1,(long double)i2,n,tol); } 00240 00241 template <class P> 00242 static bool numericallyEqual(const std::complex<P>& v1, const std::complex<P>& v2, int n, double tol=defTol<P>()) { 00243 return numericallyEqual(v1.real(), v2.real(), n, tol) 00244 && numericallyEqual(v1.imag(), v2.imag(), n, tol); 00245 } 00246 template <class P> 00247 static bool numericallyEqual(const conjugate<P>& v1, const conjugate<P>& v2, int n, double tol=defTol<P>()) { 00248 return numericallyEqual(v1.real(), v2.real(), n, tol) 00249 && numericallyEqual(v1.imag(), v2.imag(), n, tol); 00250 } 00251 template <class P> 00252 static bool numericallyEqual(const std::complex<P>& v1, const conjugate<P>& v2, int n, double tol=defTol<P>()) { 00253 return numericallyEqual(v1.real(), v2.real(), n, tol) 00254 && numericallyEqual(v1.imag(), v2.imag(), n, tol); 00255 } 00256 template <class P> 00257 static bool numericallyEqual(const conjugate<P>& v1, const std::complex<P>& v2, int n, double tol=defTol<P>()) { 00258 return numericallyEqual(v1.real(), v2.real(), n, tol) 00259 && numericallyEqual(v1.imag(), v2.imag(), n, tol); 00260 } 00261 template <class P> 00262 static bool numericallyEqual(const negator<P>& v1, const negator<P>& v2, int n, double tol=defTol<P>()) { 00263 return numericallyEqual(-v1, -v2, n, tol); // P, P 00264 } 00265 template <class P> 00266 static bool numericallyEqual(const P& v1, const negator<P>& v2, int n, double tol=defTol<P>()) { 00267 return numericallyEqual(-v1, -v2, n, tol); // P, P 00268 } 00269 template <class P> 00270 static bool numericallyEqual(const negator<P>& v1, const P& v2, int n, double tol=defTol<P>()) { 00271 return numericallyEqual(-v1, -v2, n, tol); // P, P 00272 } 00273 template <class P> 00274 static bool numericallyEqual(const negator<std::complex<P> >& v1, const conjugate<P>& v2, int n, double tol=defTol<P>()) { 00275 return numericallyEqual(-v1, -v2, n, tol); // complex, conjugate 00276 } 00277 template <class P> 00278 static bool numericallyEqual(const negator<conjugate<P> >& v1, const std::complex<P>& v2, int n, double tol=defTol<P>()) { 00279 return numericallyEqual(-v1, -v2, n, tol); // conjugate, complex 00280 } 00281 template <class P> 00282 static bool numericallyEqual(const std::complex<P>& v1, const negator<conjugate<P> >& v2, int n, double tol=defTol<P>()) { 00283 return numericallyEqual(-v1, -v2, n, tol); // complex, conjugate 00284 } 00285 template <class P> 00286 static bool numericallyEqual(const conjugate<P>& v1, const negator<std::complex<P> >& v2, int n, double tol=defTol<P>()) { 00287 return numericallyEqual(-v1, -v2, n, tol); // conjugate, complex 00288 } 00289 template <int M, class E1, int S1, class E2, int S2> 00290 static bool numericallyEqual(const Vec<M,E1,S1>& v1, const Vec<M,E2,S2>& v2, int n, double tol=(defTol2<E1,E2>())) { 00291 for (int i=0; i<M; ++i) if (!numericallyEqual(v1[i],v2[i], n, tol)) return false; 00292 return true; 00293 } 00294 template <int N, class E1, int S1, class E2, int S2> 00295 static bool numericallyEqual(const Row<N,E1,S1>& v1, const Row<N,E2,S2>& v2, int n, double tol=(defTol2<E1,E2>())) { 00296 for (int j=0; j<N; ++j) if (!numericallyEqual(v1[j],v2[j], n, tol)) return false; 00297 return true; 00298 } 00299 template <int M, int N, class E1, int CS1, int RS1, class E2, int CS2, int RS2> 00300 static bool numericallyEqual(const Mat<M,N,E1,CS1,RS1>& v1, const Mat<M,N,E2,CS2,RS2>& v2, int n, double tol=(defTol2<E1,E2>())) { 00301 for (int j=0; j<N; ++j) if (!numericallyEqual(v1(j),v2(j), n, tol)) return false; 00302 return true; 00303 } 00304 template <int N, class E1, int S1, class E2, int S2> 00305 static bool numericallyEqual(const SymMat<N,E1,S1>& v1, const SymMat<N,E2,S2>& v2, int n, double tol=(defTol2<E1,E2>())) { 00306 return numericallyEqual(v1.getAsVec(), v2.getAsVec(), n, tol); 00307 } 00308 template <class E1, class E2> 00309 static bool numericallyEqual(const VectorView_<E1>& v1, const VectorView_<E2>& v2, int n, double tol=(defTol2<E1,E2>())) { 00310 if (v1.size() != v2.size()) return false; 00311 for (int i=0; i < v1.size(); ++i) 00312 if (!numericallyEqual(v1[i], v2[i], n, tol)) return false; 00313 return true; 00314 } 00315 template <class E1, class E2> 00316 static bool numericallyEqual(const Vector_<E1>& v1, const Vector_<E2>& v2, int n, double tol=(defTol2<E1,E2>())) 00317 { return numericallyEqual((const VectorView_<E1>&)v1, (const VectorView_<E2>&)v2, n, tol); } 00318 template <class E1, class E2> 00319 static bool numericallyEqual(const Vector_<E1>& v1, const VectorView_<E2>& v2, int n, double tol=(defTol2<E1,E2>())) 00320 { return numericallyEqual((const VectorView_<E1>&)v1, (const VectorView_<E2>&)v2, n, tol); } 00321 template <class E1, class E2> 00322 static bool numericallyEqual(const VectorView_<E1>& v1, const Vector_<E2>& v2, int n, double tol=(defTol2<E1,E2>())) 00323 { return numericallyEqual((const VectorView_<E1>&)v1, (const VectorView_<E2>&)v2, n, tol); } 00324 00325 template <class E1, class E2> 00326 static bool numericallyEqual(const RowVectorView_<E1>& v1, const RowVectorView_<E2>& v2, int n, double tol=(defTol2<E1,E2>())) { 00327 if (v1.size() != v2.size()) return false; 00328 for (int i=0; i < v1.size(); ++i) 00329 if (!numericallyEqual(v1[i], v2[i], n, tol)) return false; 00330 return true; 00331 } 00332 template <class E1, class E2> 00333 static bool numericallyEqual(const RowVector_<E1>& v1, const RowVector_<E2>& v2, int n, double tol=(defTol2<E1,E2>())) 00334 { return numericallyEqual((const RowVectorView_<E1>&)v1, (const RowVectorView_<E2>&)v2, n, tol); } 00335 template <class E1, class E2> 00336 static bool numericallyEqual(const RowVector_<E1>& v1, const RowVectorView_<E2>& v2, int n, double tol=(defTol2<E1,E2>())) 00337 { return numericallyEqual((const RowVectorView_<E1>&)v1, (const RowVectorView_<E2>&)v2, n, tol); } 00338 template <class E1, class E2> 00339 static bool numericallyEqual(const RowVectorView_<E1>& v1, const RowVector_<E2>& v2, int n, double tol=(defTol2<E1,E2>())) 00340 { return numericallyEqual((const RowVectorView_<E1>&)v1, (const RowVectorView_<E2>&)v2, n, tol); } 00341 00342 template <class E1, class E2> 00343 static bool numericallyEqual(const MatrixView_<E1>& v1, const MatrixView_<E2>& v2, int n, double tol=(defTol2<E1,E2>())) { 00344 if (v1.nrow() != v2.nrow() || v1.ncol() != v2.ncol()) return false; 00345 for (int j=0; j < v1.ncol(); ++j) 00346 if (!numericallyEqual(v1(j), v2(j), n, tol)) return false; 00347 return true; 00348 } 00349 template <class E1, class E2> 00350 static bool numericallyEqual(const Matrix_<E1>& m1, const Matrix_<E2>& m2, int n, double tol=(defTol2<E1,E2>())) 00351 { return numericallyEqual((const MatrixView_<E1>&)m1, (const MatrixView_<E2>&)m2, n, tol); } 00352 template <class E1, class E2> 00353 static bool numericallyEqual(const Matrix_<E1>& m1, const MatrixView_<E2>& m2, int n, double tol=(defTol2<E1,E2>())) 00354 { return numericallyEqual((const MatrixView_<E1>&)m1, (const MatrixView_<E2>&)m2, n, tol); } 00355 template <class E1, class E2> 00356 static bool numericallyEqual(const MatrixView_<E1>& m1, const Matrix_<E2>& m2, int n, double tol=(defTol2<E1,E2>())) 00357 { return numericallyEqual((const MatrixView_<E1>&)m1, (const MatrixView_<E2>&)m2, n, tol); } 00358 00359 template <class P> 00360 static bool numericallyEqual(const Rotation_<P>& R1, const Rotation_<P>& R2, int n, double tol=defTol<P>()) { 00361 return R1.isSameRotationToWithinAngle(R2, (Real)(n*tol)); 00362 } 00363 00364 template <class P> 00365 static bool numericallyEqual(const Transform_<P>& T1, const Transform_<P>& T2, int n, double tol=defTol<P>()) { 00366 return numericallyEqual(T1.R(), T2.R(), n, tol) 00367 && numericallyEqual(T1.p(), T2.p(), n, tol); 00368 } 00369 00370 template <class P> 00371 static bool numericallyEqual(const UnitInertia_<P>& G1, const UnitInertia_<P>& G2, int n, double tol=defTol<P>()) { 00372 return numericallyEqual(G1.asSymMat33(),G2.asSymMat33(), n, tol); 00373 } 00374 00375 template <class P> 00376 static bool numericallyEqual(const Inertia_<P>& I1, const Inertia_<P>& I2, int n, double tol=defTol<P>()) { 00377 return numericallyEqual(I1.asSymMat33(),I2.asSymMat33(), n, tol); 00378 } 00379 00380 // Random numbers 00381 static Real randReal() { 00382 static Random::Uniform rand(-1,1); 00383 return rand.getValue(); 00384 } 00385 static Complex randComplex() {return Complex(randReal(),randReal());} 00386 static Conjugate randConjugate() {return Conjugate(randReal(),randReal());} 00387 static float randFloat() {return (float)randReal();} 00388 static double randDouble() {return (double)randReal();} 00389 00390 template <int M> static Vec<M> randVec() 00391 { Vec<M> v; for (int i=0; i<M; ++i) v[i]=randReal(); return v;} 00392 template <int N> static Row<N> randRow() {return ~randVec<N>();} 00393 template <int M, int N> static Mat<M,N> randMat() 00394 { Mat<M,N> m; for (int j=0; j<N; ++j) m(j)=randVec<M>(); return m;} 00395 template <int N> static SymMat<N> randSymMat() 00396 { SymMat<N> s; s.updAsVec() = randVec<N*(N+1)/2>(); return s; } 00397 00398 static Vector randVector(int m) 00399 { Vector v(m); for (int i=0; i<m; ++i) v[i]=randReal(); return v;} 00400 static Matrix randMatrix(int m, int n) 00401 { Matrix M(m,n); for (int j=0; j<n; ++j) M(j)=randVector(m); return M;} 00402 00403 static Vec3 randVec3() {return randVec<3>();} 00404 static Mat33 randMat33() {return randMat<3,3>();} 00405 static SymMat33 randSymMat33() {return randSymMat<3>();} 00406 static SpatialVec randSpatialVec() { 00407 return SpatialVec(randVec3(), randVec3()); 00408 } 00409 static SpatialMat randSpatialMat() { 00410 return SpatialMat(randMat33(), randMat33(), 00411 randMat33(), randMat33()); 00412 } 00413 static Rotation randRotation() { 00414 // Generate random angle and random axis to rotate around. 00415 return Rotation((Pi/2)*randReal(), randVec3()); 00416 } 00417 static Transform randTransform() { 00418 return Transform(randRotation(), randVec3()); 00419 } 00420 private: 00421 std::clock_t startTime; 00422 std::string testName; 00423 }; 00424 00426 class Test::Subtest { 00427 public: 00428 Subtest(const std::string& name) : subtestName(name) 00429 { 00430 char padded[128]; 00431 sprintf(padded, "%-20s", name.c_str()); 00432 paddedName = std::string(padded); 00433 std::clog << " " << paddedName << " ... " << std::flush; 00434 startTime = std::clock(); 00435 } 00436 ~Subtest() { 00437 std::clog << "done. " << paddedName << " time: " 00438 << 1000*(std::clock()-startTime)/CLOCKS_PER_SEC << "ms.\n"; 00439 } 00440 private: 00441 std::clock_t startTime; 00442 std::string subtestName; 00443 std::string paddedName; // name plus some blanks 00444 }; 00445 00446 } // namespace SimTK 00447 00449 #define SimTK_START_TEST(testName) \ 00450 SimTK::Test simtk_test_(testName); \ 00451 try { 00452 00454 #define SimTK_END_TEST() \ 00455 } catch(const std::exception& e) { \ 00456 std::cerr << "Test failed due to exception: " \ 00457 << e.what() << std::endl; \ 00458 return 1; \ 00459 } catch(...) { \ 00460 std::cerr << "Test failed due to unrecognized exception.\n"; \ 00461 return 1; \ 00462 } \ 00463 return 0; 00464 00467 #define SimTK_SUBTEST(testFunction) \ 00468 do {SimTK::Test::Subtest sub(#testFunction); (testFunction)();} while(false) 00469 00470 00471 #define SimTK_SUBTEST1(testFunction,arg1) \ 00472 do {SimTK::Test::Subtest sub(#testFunction); (testFunction)(arg1);} while(false) 00473 00474 00475 #define SimTK_SUBTEST2(testFunction,arg1,arg2) \ 00476 do {SimTK::Test::Subtest sub(#testFunction); (testFunction)(arg1,arg2);} while(false) 00477 00478 00479 #define SimTK_SUBTEST3(testFunction,arg1,arg2,arg3) \ 00480 do {SimTK::Test::Subtest sub(#testFunction); (testFunction)(arg1,arg2,arg3);} while(false) 00481 00482 00483 #define SimTK_SUBTEST4(testFunction,arg1,arg2,arg3,arg4) \ 00484 do {SimTK::Test::Subtest sub(#testFunction); (testFunction)(arg1,arg2,arg3,arg4);} while(false) 00485 00487 #define SimTK_TEST(cond) {SimTK_ASSERT_ALWAYS((cond), "Test condition failed.");} 00488 00491 #define SimTK_TEST_FAILED(msg) {SimTK_ASSERT_ALWAYS(!"Test case failed.", msg);} 00492 00496 #define SimTK_TEST_FAILED1(fmt,a1) {SimTK_ASSERT1_ALWAYS(!"Test case failed.",fmt,a1);} 00497 00501 #define SimTK_TEST_FAILED2(fmt,a1,a2) {SimTK_ASSERT2_ALWAYS(!"Test case failed.",fmt,a1,a2);} 00502 00506 #define SimTK_TEST_EQ(v1,v2) \ 00507 {SimTK_ASSERT_ALWAYS(SimTK::Test::numericallyEqual((v1),(v2),1), \ 00508 "Test values should have been numerically equivalent at default tolerance.");} 00509 00512 #define SimTK_TEST_EQ_SIZE(v1,v2,n) \ 00513 {SimTK_ASSERT1_ALWAYS(SimTK::Test::numericallyEqual((v1),(v2),(n)), \ 00514 "Test values should have been numerically equivalent at size=%d times default tolerance.",(n));} 00515 00519 #define SimTK_TEST_EQ_TOL(v1,v2,tol) \ 00520 {SimTK_ASSERT1_ALWAYS(SimTK::Test::numericallyEqual((v1),(v2),1,(tol)), \ 00521 "Test values should have been numerically equivalent at tolerance=%g.",(tol));} 00522 00526 #define SimTK_TEST_NOTEQ(v1,v2) \ 00527 {SimTK_ASSERT_ALWAYS(!SimTK::Test::numericallyEqual((v1),(v2),1), \ 00528 "Test values should NOT have been numerically equivalent (at default tolerance).");} 00529 00533 #define SimTK_TEST_NOTEQ_SIZE(v1,v2,n) \ 00534 {SimTK_ASSERT1_ALWAYS(!SimTK::Test::numericallyEqual((v1),(v2),(n)), \ 00535 "Test values should NOT have been numerically equivalent at size=%d times default tolerance.",(n));} 00536 00540 #define SimTK_TEST_NOTEQ_TOL(v1,v2,tol) \ 00541 {SimTK_ASSERT1_ALWAYS(!SimTK::Test::numericallyEqual((v1),(v2),1,(tol)), \ 00542 "Test values should NOT have been numerically equivalent at tolerance=%g.",(tol));} 00543 00545 #define SimTK_TEST_MUST_THROW(stmt) \ 00546 do {int threw=0; try {stmt;} \ 00547 catch(const std::exception&){threw=1;} \ 00548 catch(...){threw=2;} \ 00549 if (threw==0) SimTK_TEST_FAILED1("Expected statement\n----\n%s\n----\n to throw an exception but it did not.",#stmt); \ 00550 if (threw==2) SimTK_TEST_FAILED1("Expected statement\n%s\n to throw an std::exception but it threw something else.",#stmt); \ 00551 }while(false) 00552 00554 #define SimTK_TEST_MUST_THROW_EXC(stmt,exc) \ 00555 do {int threw=0; try {stmt;} \ 00556 catch(const exc&){threw=1;} \ 00557 catch(...){threw=2;} \ 00558 if (threw==0) SimTK_TEST_FAILED1("Expected statement\n----\n%s\n----\n to throw an exception but it did not.",#stmt); \ 00559 if (threw==2) SimTK_TEST_FAILED2("Expected statement\n----\n%s\n----\n to throw exception type %s but it threw something else.",#stmt,#exc); \ 00560 }while(false) 00561 00563 #define SimTK_TEST_MAY_THROW(stmt) \ 00564 do {int threw=0; try {stmt;} \ 00565 catch(const std::exception&){threw=1;} \ 00566 catch(...){threw=2;} \ 00567 if (threw==2) SimTK_TEST_FAILED1("Expected statement\n%s\n to throw an std::exception but it threw something else.",#stmt); \ 00568 }while(false) 00569 00571 #define SimTK_TEST_MAY_THROW_EXC(stmt,exc) \ 00572 do {int threw=0; try {stmt;} \ 00573 catch(const exc&){threw=1;} \ 00574 catch(...){threw=2;} \ 00575 if (threw==2) SimTK_TEST_FAILED2("Expected statement\n----\n%s\n----\n to throw exception type %s but it threw something else.",#stmt,#exc); \ 00576 }while(false) 00577 00578 // When we're only required to throw in Debug, we have to suppress the 00579 // test case altogether in Release because it may cause damage. 00580 #if defined(NDEBUG) 00581 00582 00583 #define SimTK_TEST_MUST_THROW_DEBUG(stmt) 00584 00585 00586 #define SimTK_TEST_MUST_THROW_EXC_DEBUG(stmt,exc) 00587 #else 00588 00589 00590 #define SimTK_TEST_MUST_THROW_DEBUG(stmt) SimTK_TEST_MUST_THROW(stmt) 00591 00592 00593 #define SimTK_TEST_MUST_THROW_EXC_DEBUG(stmt,exc) \ 00594 SimTK_TEST_MUST_THROW_EXC(stmt,exc) 00595 #endif 00596 00597 00598 00599 00600 // End of Regression testing group. 00602 00603 #endif // SimTK_SimTKCOMMON_TESTING_H_