00001
00002 #include "Math/SMatrix.h"
00003 #include "TMatrixTSym.h"
00004 #include "TDecompChol.h"
00005 #include "TDecompBK.h"
00006
00007 #include "TRandom.h"
00008 #include <vector>
00009 #include <string>
00010 #include <iostream>
00011 #include <cmath>
00012 #include <limits>
00013
00014 #include "TStopwatch.h"
00015
00016
00017 #ifndef N
00018 #define N 5
00019 #endif
00020
00021 bool doSelfTest = true;
00022
00023
00024 namespace test {
00025
00026
00027 #ifdef REPORT_TIME
00028 void reportTime( std::string s, double time);
00029 #endif
00030
00031 void printTime(TStopwatch & time, std::string s) {
00032 int pr = std::cout.precision(8);
00033 std::cout << s << "\t" << " time = " << time.RealTime() << "\t(sec)\t"
00034
00035 << std::endl;
00036 std::cout.precision(pr);
00037 }
00038
00039
00040
00041 class Timer {
00042
00043 public:
00044
00045 Timer(const std::string & s = "") : fName(s), fTime(0)
00046 {
00047 fWatch.Start();
00048 }
00049 Timer(double & t, const std::string & s = "") : fName(s), fTime(&t)
00050 {
00051 fWatch.Start();
00052 }
00053
00054 ~Timer() {
00055 fWatch.Stop();
00056 printTime(fWatch,fName);
00057 #ifdef REPORT_TIME
00058
00059 reportTime(fName, fWatch.RealTime() );
00060 #endif
00061 if (fTime) *fTime += fWatch.RealTime();
00062 }
00063
00064
00065 private:
00066
00067 std::string fName;
00068 double * fTime;
00069 TStopwatch fWatch;
00070
00071 };
00072 }
00073
00074 using namespace ROOT::Math;
00075
00076
00077
00078 typedef SMatrix<double,N,N, MatRepSym<double,N> > SymMatrix;
00079
00080
00081 template<class M>
00082 M * createMatrix() {
00083 return new M();
00084 }
00085
00086 template<>
00087 TMatrixTSym<double> * createMatrix<TMatrixTSym<double> >() {
00088 return new TMatrixTSym<double>(N);
00089 }
00090
00091
00092 template<class M>
00093 void printMatrix(const M & m) {
00094 std::cout << m << std::endl;
00095 }
00096 template<>
00097 void printMatrix<TMatrixTSym<double> >(const TMatrixTSym<double> & m ) {
00098 m.Print();
00099 }
00100
00101
00102 template<class M>
00103 void genMatrix(M & m ) {
00104 TRandom & r = *gRandom;
00105
00106 for (int i = 0; i < N; ++i) {
00107 double maxVal = i*10000/(N-1) + 1;
00108 m(i,i) = r.Uniform(0, maxVal);
00109 }
00110 for (int i = 0; i < N; ++i) {
00111 for (int j = 0; j < i; ++j) {
00112 double v = 0.3*std::sqrt( m(i,i) * m(j,j) );
00113 m(i,j) = r.Uniform(0, v);
00114 m(j,i) = m(i,j);
00115 }
00116 }
00117 }
00118
00119
00120 template<class M>
00121 void generate(std::vector<M*> & v) {
00122 int n = v.size();
00123 gRandom->SetSeed(111);
00124 for (int i = 0; i < n; ++i) {
00125 v[i] = createMatrix<M>();
00126 genMatrix(*v[i] );
00127 }
00128 }
00129
00130
00131 struct Choleski {};
00132 struct BK {};
00133 struct QR {};
00134 struct Cramer {};
00135 struct Default {};
00136
00137 template <class M, class Type>
00138 struct TestInverter {
00139 static bool Inv ( const M & , M & ) { return false;}
00140 static bool Inv2 ( M & ) { return false;}
00141 };
00142
00143 template <>
00144 struct TestInverter<SymMatrix, Choleski> {
00145 static bool Inv ( const SymMatrix & m, SymMatrix & result ) {
00146 int ifail = 0;
00147 result = m.InverseChol(ifail);
00148 return ifail == 0;
00149 }
00150 static bool Inv2 ( SymMatrix & m ) {
00151 return m.InvertChol();
00152 }
00153 };
00154
00155 template <>
00156 struct TestInverter<SymMatrix, BK> {
00157 static bool Inv ( const SymMatrix & m, SymMatrix & result ) {
00158 int ifail = 0;
00159 result = m.Inverse(ifail);
00160 return ifail==0;
00161 }
00162 static bool Inv2 ( SymMatrix & m ) {
00163 return m.Invert();
00164 }
00165 };
00166
00167 template <>
00168 struct TestInverter<SymMatrix, Cramer> {
00169 static bool Inv ( const SymMatrix & m, SymMatrix & result ) {
00170 int ifail = 0;
00171 result = m.InverseFast(ifail);
00172 return ifail==0;
00173 }
00174 static bool Inv2 ( SymMatrix & m ) {
00175 return m.InvertFast();
00176 }
00177 };
00178
00179 #ifdef LATER
00180 template <>
00181 struct TestInverter<SymMatrix, QR> {
00182 static bool Inv ( const SymMatrix & m, SymMatrix & result ) {
00183 ROOT::Math::QRDecomposition<double> d;
00184 int ifail = 0;
00185 result = m.InverseFast(ifail);
00186 return ifail==0;
00187 }
00188 };
00189 #endif
00190
00191
00192
00193 template <>
00194 struct TestInverter<TMatrixDSym, Default> {
00195 static bool Inv ( const TMatrixDSym & m, TMatrixDSym & result ) {
00196 result = m;
00197 result.Invert();
00198 return true;
00199 }
00200 static bool Inv2 ( TMatrixDSym & m ) {
00201 m.Invert();
00202 return true;
00203 }
00204 };
00205
00206 template <>
00207 struct TestInverter<TMatrixDSym, Cramer> {
00208 static bool Inv ( const TMatrixDSym & m, TMatrixDSym & result ) {
00209 result = m;
00210 result.InvertFast();
00211 return true;
00212 }
00213 static bool Inv2 ( TMatrixDSym & m ) {
00214 m.InvertFast();
00215 return true;
00216 }
00217 };
00218
00219 template <>
00220 struct TestInverter<TMatrixDSym, Choleski> {
00221 static bool Inv ( const TMatrixDSym & m, TMatrixDSym & result ) {
00222 TDecompChol chol(m);
00223 if (!chol.Decompose() ) return false;
00224 chol.Invert(result);
00225 return true;
00226 }
00227 };
00228
00229 template <>
00230 struct TestInverter<TMatrixDSym, BK> {
00231 static bool Inv ( const TMatrixDSym & m, TMatrixDSym & result ) {
00232 TDecompBK d(m);
00233 if (!d.Decompose() ) return false;
00234 d.Invert(result);
00235 return true;
00236 }
00237 };
00238
00239
00240 template<class M, class T>
00241 double invert( const std::vector<M* > & matlist, double & time,std::string s) {
00242 M result = *(matlist.front());
00243 test::Timer t(time,s);
00244 int nloop = matlist.size();
00245 double sum = 0;
00246 for (int l = 0; l < nloop; l++)
00247 {
00248 const M & m = *(matlist[l]);
00249 bool ok = TestInverter<M,T>::Inv(m,result);
00250 if (!ok) {
00251 std::cout << "inv failed for matrix " << l << std::endl;
00252 printMatrix<M>( m);
00253 return -1;
00254 }
00255 sum += result(0,1);
00256 }
00257 return sum;
00258 }
00259
00260
00261
00262 template<class M, class T>
00263 double invert2( const std::vector<M* > & matlist, double & time,std::string s) {
00264
00265
00266 int nloop = matlist.size();
00267 std::vector<M *> vmat(nloop);
00268 for (int l = 0; l < nloop; l++)
00269 {
00270 vmat[l] = new M( *matlist[l] );
00271 }
00272
00273 test::Timer t(time,s);
00274 double sum = 0;
00275 for (int l = 0; l < nloop; l++)
00276 {
00277 M & m = *(vmat[l]);
00278 bool ok = TestInverter<M,T>::Inv2(m);
00279 if (!ok) {
00280 std::cout << "inv failed for matrix " << l << std::endl;
00281 printMatrix<M>( m);
00282 return -1;
00283 }
00284 sum += m(0,1);
00285 }
00286 return sum;
00287 }
00288
00289 bool equal(double d1, double d2, double stol = 10000) {
00290 std::cout.precision(12);
00291 double eps = stol * std::numeric_limits<double>::epsilon();
00292 if ( std::abs(d1) > 0 && std::abs(d2) > 0 )
00293 return ( std::abs( d1-d2) < eps * std::max(std::abs(d1), std::abs(d2) ) );
00294 else if ( d1 == 0 )
00295 return std::abs(d2) < eps;
00296 else
00297 return std::abs(d1) < eps;
00298 }
00299
00300
00301 bool stressSymPosInversion(int n, bool selftest ) {
00302
00303
00304
00305 std::vector<SymMatrix *> v1(n);
00306 generate(v1);
00307 std::vector<TMatrixDSym *> v2(n);
00308 generate(v2);
00309
00310
00311 bool iret = true;
00312 double time = 0;
00313 double s1 = invert<SymMatrix, Choleski> (v1, time,"SMatrix Chol");
00314 double s2 = invert<SymMatrix, BK> (v1, time,"SMatrix BK");
00315 double s3 = invert<SymMatrix, Cramer> (v1, time,"SMatrix Cram");
00316 bool ok = ( equal(s1,s2) && equal(s1,s3) );
00317 if (!ok) {
00318 std::cout << "result SMatrix choleski " << s1 << " BK " << s2 << " cramer " << s3 << std::endl;
00319 std::cerr <<"Error: inversion test for SMatrix FAILED ! " << std::endl;
00320 }
00321 iret &= ok;
00322 std::cout << std::endl;
00323
00324 double m1 = invert<TMatrixDSym, Choleski> (v2, time,"TMatrix Chol");
00325 double m2 = invert<TMatrixDSym, BK> (v2, time,"TMatrix BK");
00326 double m3 = invert<TMatrixDSym, Cramer> (v2, time,"TMatrix Cram");
00327 double m4 = invert<TMatrixDSym, Default> (v2, time,"TMatrix Def");
00328
00329 ok = ( equal(m1,m2) && equal(m1,m3) && equal(m1,m4) );
00330 if (!ok) {
00331 std::cout << "result TMatrix choleski " << m1 << " BK " << m2
00332 << " cramer " << m3 << " default " << m4 << std::endl;
00333 std::cerr <<"Error: inversion test for TMatrix FAILED ! " << std::endl;
00334 }
00335 iret &= ok;
00336 std::cout << std::endl;
00337
00338
00339
00340 if (selftest) {
00341 std::cout << "\n - self inversion test \n";
00342 double s11 = invert2<SymMatrix, Choleski> (v1, time,"SMatrix Chol");
00343 double s12 = invert2<SymMatrix, BK> (v1, time,"SMatrix BK");
00344 double s13 = invert2<SymMatrix, Cramer> (v1, time,"SMatrix Cram");
00345 ok = ( equal(s11,s12) && equal(s11,s13) );
00346 if (!ok) {
00347 std::cout << "result SMatrix choleski " << s11 << " BK " << s12 << " cramer " << s13 << std::endl;
00348 std::cerr <<"Error: self inversion test for SMatrix FAILED ! " << std::endl;
00349 }
00350 iret &= ok;
00351 std::cout << std::endl;
00352
00353 double m13 = invert2<TMatrixDSym, Cramer> (v2, time,"TMatrix Cram");
00354 double m14 = invert2<TMatrixDSym, Default> (v2, time,"TMatrix Def");
00355 ok = ( equal(m13,m14) );
00356 if (!ok) {
00357 std::cout << "result TMatrix cramer " << m13 << " default " << m14 << std::endl;
00358 std::cerr <<"Error: self inversion test for TMatrix FAILED ! " << std::endl;
00359 }
00360 iret &= ok;
00361 std::cout << std::endl;
00362 }
00363
00364 return iret;
00365 }
00366
00367 int testInversion(int n = 100000) {
00368 bool ok = stressSymPosInversion(n, doSelfTest);
00369 std::cerr << "Test inversion of positive defined matrix ....... ";
00370 if (ok) std::cerr << "OK \n";
00371 else std::cerr << "FAILED \n";
00372 return (ok) ? 0 : -1;
00373 }
00374
00375 int main() {
00376 return testInversion();
00377 }