Tools.cxx

Go to the documentation of this file.
00001 // @(#)root/tmva $Id: Tools.cxx 37986 2011-02-04 21:42:15Z pcanal $
00002 // Author: Andreas Hoecker, Peter Speckmayer, Joerg Stelzer, Helge Voss
00003 
00004 /**********************************************************************************
00005  * Project: TMVA - a Root-integrated toolkit for multivariate data analysis       *
00006  * Package: TMVA                                                                  *
00007  * Class  : Tools                                                                 *
00008  * Web    : http://tmva.sourceforge.net                                           *
00009  *                                                                                *
00010  * Description:                                                                   *
00011  *      Implementation (see header for description)                               *
00012  *                                                                                *
00013  * Authors (alphabetical):                                                        *
00014  *      Andreas Hoecker <Andreas.Hocker@cern.ch> - CERN, Switzerland              *
00015  *      Peter Speckmayer <Peter.Speckmayer@cern.ch> - CERN, Switzerland           *
00016  *      Helge Voss      <Helge.Voss@cern.ch>     - MPI-K Heidelberg, Germany      *
00017  *      Kai Voss        <Kai.Voss@cern.ch>       - U. of Victoria, Canada         *
00018  *                                                                                *
00019  * Copyright (c) 2005:                                                            *
00020  *      CERN, Switzerland                                                         *
00021  *      U. of Victoria, Canada                                                    *
00022  *      MPI-K Heidelberg, Germany                                                 *
00023  *                                                                                *
00024  * Redistribution and use in source and binary forms, with or without             *
00025  * modification, are permitted according to the terms listed in LICENSE           *
00026  * (http://ttmva.sourceforge.net/LICENSE)                                         *
00027  **********************************************************************************/
00028 
00029 #include <algorithm>
00030 #include <cstdlib>
00031 
00032 #include "TObjString.h"
00033 #include "TMath.h"
00034 #include "TString.h"
00035 #include "TTree.h"
00036 #include "TLeaf.h"
00037 #include "TH1.h"
00038 #include "TH2.h"
00039 #include "TList.h"
00040 #include "TSpline.h"
00041 #include "TVector.h"
00042 #include "TMatrixD.h"
00043 #include "TMatrixDSymEigen.h"
00044 #include "TVectorD.h"
00045 #include "TTreeFormula.h"
00046 #include "TXMLEngine.h"
00047 #include "TROOT.h"
00048 #include "TMatrixDSymEigen.h"
00049 
00050 #ifndef ROOT_TMVA_Tools
00051 #include "TMVA/Tools.h"
00052 #endif
00053 #ifndef ROOT_TMVA_Config
00054 #include "TMVA/Config.h"
00055 #endif
00056 #ifndef ROOT_TMVA_Event
00057 #include "TMVA/Event.h"
00058 #endif
00059 #ifndef ROOT_TMVA_Version
00060 #include "TMVA/Version.h"
00061 #endif
00062 #ifndef ROOT_TMVA_PDF
00063 #include "TMVA/PDF.h"
00064 #endif
00065 #ifndef ROOT_TMVA_MsgLogger
00066 #include "TMVA/MsgLogger.h"
00067 #endif
00068 
00069 using namespace std;
00070 
00071 TMVA::Tools* TMVA::Tools::fgTools = 0;
00072 TMVA::Tools& TMVA::gTools()                 { return TMVA::Tools::Instance(); }
00073 TMVA::Tools& TMVA::Tools::Instance()        { return fgTools?*(fgTools): *(fgTools = new Tools()); }
00074 void         TMVA::Tools::DestroyInstance() { if (fgTools != 0) { delete fgTools; fgTools=0; } }
00075 
00076 //_______________________________________________________________________
00077 TMVA::Tools::Tools() :
00078    fRegexp("$&|!%^&()'<>?= "),
00079    fLogger(new MsgLogger("Tools")),
00080    fXMLEngine(new TXMLEngine())
00081 {
00082    // constructor
00083 }
00084 
00085 //_______________________________________________________________________
00086 TMVA::Tools::~Tools()
00087 {
00088    // destructor
00089    delete fLogger;
00090    delete fXMLEngine;
00091 }
00092 
00093 //_______________________________________________________________________
00094 Double_t TMVA::Tools::NormVariable( Double_t x, Double_t xmin, Double_t xmax )
00095 {
00096    // normalise to output range: [-1, 1]
00097    return 2*(x - xmin)/(xmax - xmin) - 1.0;
00098 }
00099 
00100 //_______________________________________________________________________
00101 Double_t TMVA::Tools::GetSeparation( TH1* S, TH1* B ) const
00102 {
00103    // compute "separation" defined as
00104    // <s2> = (1/2) Int_-oo..+oo { (S^2(x) - B^2(x))/(S(x) + B(x)) dx }
00105    Double_t separation = 0;
00106 
00107    // sanity checks
00108    // signal and background histograms must have same number of bins and 
00109    // same limits
00110    if ((S->GetNbinsX() != B->GetNbinsX()) || (S->GetNbinsX() <= 0)) {
00111       Log() << kFATAL << "<GetSeparation> signal and background"
00112             << " histograms have different number of bins: " 
00113             << S->GetNbinsX() << " : " << B->GetNbinsX() << Endl;
00114    }
00115 
00116    if (S->GetXaxis()->GetXmin() != B->GetXaxis()->GetXmin() || 
00117        S->GetXaxis()->GetXmax() != B->GetXaxis()->GetXmax() || 
00118        S->GetXaxis()->GetXmax() <= S->GetXaxis()->GetXmin()) {
00119       Log() << kINFO << S->GetXaxis()->GetXmin() << " " << B->GetXaxis()->GetXmin() 
00120             << " " << S->GetXaxis()->GetXmax() << " " << B->GetXaxis()->GetXmax() 
00121             << " " << S->GetXaxis()->GetXmax() << " " << S->GetXaxis()->GetXmin() << Endl;
00122       Log() << kFATAL << "<GetSeparation> signal and background"
00123             << " histograms have different or invalid dimensions:" << Endl;
00124    }
00125 
00126    Int_t    nstep  = S->GetNbinsX();
00127    Double_t intBin = (S->GetXaxis()->GetXmax() - S->GetXaxis()->GetXmin())/nstep;
00128    Double_t nS     = S->GetSumOfWeights()*intBin;
00129    Double_t nB     = B->GetSumOfWeights()*intBin;
00130 
00131    if (nS > 0 && nB > 0) {
00132       for (Int_t bin=0; bin<nstep; bin++) {
00133          Double_t s = S->GetBinContent( bin )/Double_t(nS);
00134          Double_t b = B->GetBinContent( bin )/Double_t(nB);
00135          // separation
00136          if (s + b > 0) separation += 0.5*(s - b)*(s - b)/(s + b);
00137       }
00138       separation *= intBin;
00139    }
00140    else {
00141       Log() << kWARNING << "<GetSeparation> histograms with zero entries: " 
00142             << nS << " : " << nB << " cannot compute separation"
00143             << Endl;
00144       separation = 0;
00145    }
00146 
00147    return separation;
00148 }
00149 
00150 //_______________________________________________________________________
00151 Double_t TMVA::Tools::GetSeparation( const PDF& pdfS, const PDF& pdfB ) const
00152 {
00153    // compute "separation" defined as
00154    // <s2> = (1/2) Int_-oo..+oo { (S(x)2 - B(x)2)/(S(x) + B(x)) dx }
00155 
00156    Double_t xmin = pdfS.GetXmin();
00157    Double_t xmax = pdfS.GetXmax();
00158    // sanity check
00159    if (xmin != pdfB.GetXmin() || xmax != pdfB.GetXmax()) {
00160       Log() << kFATAL << "<GetSeparation> Mismatch in PDF limits: "
00161             << xmin << " " << pdfB.GetXmin() << xmax << " " << pdfB.GetXmax()  << Endl;
00162    }
00163 
00164    Double_t separation = 0;
00165    Int_t    nstep      = 100;
00166    Double_t intBin     = (xmax - xmin)/Double_t(nstep);
00167    for (Int_t bin=0; bin<nstep; bin++) {
00168       Double_t x = (bin + 0.5)*intBin + xmin;
00169       Double_t s = pdfS.GetVal( x );
00170       Double_t b = pdfB.GetVal( x );
00171       // separation
00172       if (s + b > 0) separation += (s - b)*(s - b)/(s + b);
00173    }
00174    separation *= (0.5*intBin);
00175 
00176    return separation;
00177 }
00178 
00179 //_______________________________________________________________________
00180 void TMVA::Tools::ComputeStat( const std::vector<TMVA::Event*>& events, std::vector<Float_t>* valVec,
00181                                Double_t& meanS, Double_t& meanB,
00182                                Double_t& rmsS,  Double_t& rmsB,
00183                                Double_t& xmin,  Double_t& xmax,
00184                                Int_t signalClass, Bool_t  norm )
00185 {
00186    // sanity check
00187    if (0 == valVec) 
00188       Log() << kFATAL << "<Tools::ComputeStat> value vector is zero pointer" << Endl;
00189    
00190    if ( events.size() != valVec->size() ) 
00191       Log() << kFATAL << "<Tools::ComputeStat> event and value vector have different lengths " 
00192             << events.size() << "!=" << valVec->size() << Endl;
00193 
00194    Long64_t entries = valVec->size();
00195 
00196    // first fill signal and background in arrays before analysis
00197    Double_t* varVecS  = new Double_t[entries];
00198    Double_t* varVecB  = new Double_t[entries];
00199    xmin               = +DBL_MAX;
00200    xmax               = -DBL_MAX;
00201    Long64_t nEventsS  = -1;
00202    Long64_t nEventsB  = -1;
00203    Double_t xmin_ = 0, xmax_ = 0;
00204 
00205    if (norm) {
00206       xmin_ = *std::min( valVec->begin(), valVec->end() );
00207       xmax_ = *std::max( valVec->begin(), valVec->end() );
00208    }
00209 
00210    for (Int_t ievt=0; ievt<entries; ievt++) {
00211       Double_t theVar = (*valVec)[ievt];
00212       if (norm) theVar = Tools::NormVariable( theVar, xmin_, xmax_ );
00213 
00214       if (Int_t(events[ievt]->GetClass()) == signalClass ){
00215          varVecS[++nEventsS] = theVar; // this is signal
00216       }
00217       else {
00218          varVecB[++nEventsB] = theVar; // this is background
00219       }
00220 
00221       if (theVar > xmax) xmax = theVar;
00222       if (theVar < xmin) xmin = theVar;
00223    }
00224    ++nEventsS;
00225    ++nEventsB;
00226 
00227    // basic statistics
00228    meanS = TMath::Mean( nEventsS, varVecS );
00229    meanB = TMath::Mean( nEventsB, varVecB );
00230    rmsS  = TMath::RMS ( nEventsS, varVecS );
00231    rmsB  = TMath::RMS ( nEventsB, varVecB );
00232 
00233    delete [] varVecS;
00234    delete [] varVecB;
00235 }
00236 
00237 //_______________________________________________________________________
00238 TMatrixD* TMVA::Tools::GetSQRootMatrix( TMatrixDSym* symMat )
00239 {
00240    // square-root of symmetric matrix
00241    // of course the resulting sqrtMat is also symmetric, but it's easier to
00242    // treat it as a general matrix
00243    Int_t n = symMat->GetNrows();
00244 
00245    // compute eigenvectors
00246    TMatrixDSymEigen* eigen = new TMatrixDSymEigen( *symMat );
00247 
00248    // D = ST C S
00249    TMatrixD* si = new TMatrixD( eigen->GetEigenVectors() );
00250    TMatrixD* s  = new TMatrixD( *si ); // copy
00251    si->Transpose( *si ); // invert (= transpose)
00252 
00253    // diagonal matrices
00254    TMatrixD* d = new TMatrixD( n, n);
00255    d->Mult( (*si), (*symMat) ); (*d) *= (*s);
00256 
00257    // sanity check: matrix must be diagonal and positive definit
00258    Int_t i, j;
00259    Double_t epsilon = 1.0e-8;
00260    for (i=0; i<n; i++) {
00261       for (j=0; j<n; j++) {
00262          if ((i != j && TMath::Abs((*d)(i,j))/TMath::Sqrt((*d)(i,i)*(*d)(j,j)) > epsilon) ||
00263              (i == j && (*d)(i,i) < 0)) {
00264             //d->Print();
00265             Log() << kWARNING << "<GetSQRootMatrix> error in matrix diagonalization; printed S and B" << Endl;
00266          }
00267       }
00268    }
00269 
00270    // make exactly diagonal
00271    for (i=0; i<n; i++) for (j=0; j<n; j++) if (j != i) (*d)(i,j) = 0;
00272 
00273    // compute the square-root C' of covariance matrix: C = C'*C'
00274    for (i=0; i<n; i++) (*d)(i,i) = TMath::Sqrt((*d)(i,i));
00275 
00276    TMatrixD* sqrtMat = new TMatrixD( n, n );
00277    sqrtMat->Mult( (*s), (*d) );
00278    (*sqrtMat) *= (*si);
00279 
00280    // invert square-root matrices
00281    sqrtMat->Invert();
00282 
00283    delete eigen;
00284    delete s;
00285    delete si;
00286    delete d;
00287 
00288    return sqrtMat;
00289 }
00290 
00291 //_______________________________________________________________________
00292 const TMatrixD* TMVA::Tools::GetCorrelationMatrix( const TMatrixD* covMat )
00293 {
00294    // turns covariance into correlation matrix   
00295    if (covMat == 0) return 0;
00296 
00297    // sanity check
00298    Int_t nvar = covMat->GetNrows();
00299    if (nvar != covMat->GetNcols()) 
00300       Log() << kFATAL << "<GetCorrelationMatrix> input matrix not quadratic" << Endl;
00301    
00302    TMatrixD* corrMat = new TMatrixD( nvar, nvar );
00303 
00304    for (Int_t ivar=0; ivar<nvar; ivar++) {
00305       for (Int_t jvar=0; jvar<nvar; jvar++) {
00306          if (ivar != jvar) {
00307             Double_t d = (*covMat)(ivar, ivar)*(*covMat)(jvar, jvar);
00308             if (d > 1E-20) (*corrMat)(ivar, jvar) = (*covMat)(ivar, jvar)/TMath::Sqrt(d);
00309             else {
00310                Log() << kWARNING << "<GetCorrelationMatrix> zero variances for variables "
00311                      << "(" << ivar << ", " << jvar << ")" << Endl;
00312                (*corrMat)(ivar, jvar) = 0;
00313             }
00314             if (TMath::Abs( (*corrMat)(ivar,jvar))  > 1){
00315                Log() << kWARNING
00316                      <<  " Element  corr("<<ivar<<","<<ivar<<")=" << (*corrMat)(ivar,jvar)  
00317                      << " sigma2="<<d
00318                      << " cov("<<ivar<<","<<ivar<<")=" <<(*covMat)(ivar, ivar)
00319                      << " cov("<<jvar<<","<<jvar<<")=" <<(*covMat)(jvar, jvar)
00320                      << Endl; 
00321                
00322             }
00323          }
00324          else (*corrMat)(ivar, ivar) = 1.0;
00325       }
00326    }
00327 
00328    return corrMat;
00329 }
00330 
00331 //_______________________________________________________________________
00332 TH1* TMVA::Tools::projNormTH1F( TTree* theTree, const TString& theVarName,
00333                                 const TString& name, Int_t nbins,
00334                                 Double_t xmin, Double_t xmax, const TString& cut )
00335 {
00336    // projects variable from tree into normalised histogram
00337  
00338    // needed because of ROOT bug (feature) that excludes events that have value == xmax
00339    xmax += 0.00001; 
00340    
00341    TH1* hist = new TH1F( name, name, nbins, xmin, xmax );
00342    hist->Sumw2(); // enable quadratic errors
00343    theTree->Project( name, theVarName, cut );
00344    NormHist( hist );
00345    return hist;
00346 }
00347 
00348 //_______________________________________________________________________
00349 Double_t TMVA::Tools::NormHist( TH1* theHist, Double_t norm )
00350 {
00351    // normalises histogram
00352    if (!theHist) return 0;
00353 
00354    if (theHist->GetSumw2N() == 0) theHist->Sumw2();
00355    if (theHist->GetSumOfWeights() != 0) {
00356       Double_t   w  = ( theHist->GetSumOfWeights()
00357                         *(theHist->GetXaxis()->GetXmax() - theHist->GetXaxis()->GetXmin())/theHist->GetNbinsX() );
00358       if (w > 0) theHist->Scale( norm/w );
00359       return w;
00360    }
00361 
00362    return 1.0;
00363 }
00364 
00365 //_______________________________________________________________________
00366 TList* TMVA::Tools::ParseFormatLine( TString formatString, const char* sep )
00367 {
00368    // Parse the string and cut into labels separated by ":"
00369    TList*   labelList = new TList();
00370    labelList->SetOwner();
00371    while (formatString.First(sep)==0) formatString.Remove(0,1); // remove initial separators
00372 
00373    while (formatString.Length()>0) {
00374       if (formatString.First(sep) == -1) { // no more separator
00375          labelList->Add(new TObjString(formatString.Data()));
00376          formatString="";
00377          break;
00378       }
00379 
00380       Ssiz_t posSep = formatString.First(sep);
00381       labelList->Add(new TObjString(TString(formatString(0,posSep)).Data()));
00382       formatString.Remove(0,posSep+1);
00383       
00384       while (formatString.First(sep)==0) formatString.Remove(0,1); // remove additional separators
00385       
00386    }
00387    return labelList;                                                 
00388 }
00389 
00390 //_______________________________________________________________________
00391 vector<Int_t>* TMVA::Tools::ParseANNOptionString( TString theOptions, Int_t nvar,
00392                                                   vector<Int_t>* nodes )
00393 {
00394    // parse option string for ANN methods
00395    // default settings (should be defined in theOption string)
00396    TList* list  = TMVA::Tools::ParseFormatLine( theOptions, ":" );
00397 
00398    // format and syntax of option string: "3000:N:N+2:N-3:6"
00399    //
00400    // where:
00401    //        3000 - number of training cycles (epochs)
00402    //        N    - number of nodes in first hidden layer, where N is the number
00403    //               of discriminating variables used (note that the first ANN
00404    //               layer necessarily has N nodes, and hence is not given).
00405    //        N+2  - number of nodes in 2nd hidden layer (2 nodes more than
00406    //               number of variables)
00407    //        N-3  - number of nodes in 3rd hidden layer (3 nodes less than
00408    //               number of variables)
00409    //        6    - 6 nodes in last (4th) hidden layer (note that the last ANN
00410    //               layer in MVA has 2 nodes, each one for signal and background
00411    //               classes)
00412 
00413    // sanity check
00414    if (list->GetSize() < 1) {
00415       Log() << kFATAL << "<ParseANNOptionString> unrecognized option string: " << theOptions << Endl;
00416    }
00417 
00418    // add number of cycles
00419    nodes->push_back( atoi( ((TObjString*)list->At(0))->GetString() ) );
00420 
00421    Int_t a;
00422    if (list->GetSize() > 1) {
00423       for (Int_t i=1; i<list->GetSize(); i++) {
00424          TString s = ((TObjString*)list->At(i))->GetString();
00425          s.ToUpper();
00426          if (s(0) == 'N')  {
00427             if (s.Length() > 1) nodes->push_back( nvar + atoi(&s[1]) );
00428             else                nodes->push_back( nvar );
00429          }
00430          else if ((a = atoi( s )) > 0) nodes->push_back( atoi(s ) );
00431          else {
00432             Log() << kFATAL << "<ParseANNOptionString> unrecognized option string: " << theOptions << Endl;
00433          }
00434       }
00435    }
00436 
00437    return nodes;
00438 }
00439 
00440 Bool_t TMVA::Tools::CheckSplines( const TH1* theHist, const TSpline* theSpline )
00441 {
00442    // check quality of splining by comparing splines and histograms in each bin
00443    const Double_t sanityCrit = 0.01; // relative deviation
00444 
00445    Bool_t retval = kTRUE;
00446    for (Int_t ibin=1; ibin<=theHist->GetNbinsX(); ibin++) {
00447       Double_t x  = theHist->GetBinCenter( ibin );
00448       Double_t yh = theHist->GetBinContent( ibin ); // the histogram output
00449       Double_t ys = theSpline->Eval( x );           // the spline output
00450 
00451       if (ys + yh > 0) {
00452          Double_t dev = 0.5*(ys - yh)/(ys + yh);
00453          if (TMath::Abs(dev) > sanityCrit) {
00454             Log() << kFATAL << "<CheckSplines> Spline failed sanity criterion; "
00455                   << " relative deviation from histogram: " << dev
00456                   << " in (bin, value): (" << ibin << ", " << x << ")" << Endl;
00457             retval = kFALSE;
00458          }
00459       }
00460    }
00461 
00462    return retval;
00463 }
00464 
00465 //_______________________________________________________________________
00466 std::vector<Double_t> TMVA::Tools::MVADiff( std::vector<Double_t>& a, std::vector<Double_t>& b )
00467 {
00468    // computes difference between two vectors
00469    if (a.size() != b.size()) {
00470       throw;
00471    }
00472    vector<Double_t> result(a.size());
00473    for (UInt_t i=0; i<a.size();i++) result[i]=a[i]-b[i];
00474    return result;
00475 }
00476 
00477 //_______________________________________________________________________
00478 void TMVA::Tools::Scale( std::vector<Double_t>& v, Double_t f )
00479 {
00480    // scales double vector
00481    for (UInt_t i=0; i<v.size();i++) v[i]*=f;
00482 }
00483 
00484 //_______________________________________________________________________
00485 void TMVA::Tools::Scale( std::vector<Float_t>& v, Float_t f )
00486 {
00487    // scales float vector
00488    for (UInt_t i=0; i<v.size();i++) v[i]*=f;
00489 }
00490 
00491 //_______________________________________________________________________
00492 void TMVA::Tools::UsefulSortAscending( std::vector<vector<Double_t> >& v, std::vector<TString>* vs ){
00493    // sort 2D vector (AND in parallel a TString vector) in such a way 
00494    // that the "first vector is sorted" and the other vectors are reshuffled
00495    // in the same way as necessary to have the first vector sorted. 
00496    // I.e. the correlation between the elements is kept.
00497    UInt_t nArrays=v.size();
00498    Double_t temp;
00499    if (nArrays > 0) {
00500       UInt_t sizeofarray=v[0].size();
00501       for (UInt_t i=0; i<sizeofarray; i++) {
00502          for (UInt_t j=sizeofarray-1; j>i; j--) {
00503             if (v[0][j-1] > v[0][j]) {
00504                for (UInt_t k=0; k< nArrays; k++) {
00505                   temp = v[k][j-1]; v[k][j-1] = v[k][j]; v[k][j] = temp;
00506                }
00507                if (NULL != vs) {
00508                   TString temps = (*vs)[j-1]; (*vs)[j-1] = (*vs)[j]; (*vs)[j] = temps;
00509                }
00510             }
00511          }
00512       }
00513    }
00514 }
00515 
00516 //_______________________________________________________________________
00517 void TMVA::Tools::UsefulSortDescending( std::vector<std::vector<Double_t> >& v, std::vector<TString>* vs )
00518 {
00519    // sort 2D vector (AND in parallel a TString vector) in such a way 
00520    // that the "first vector is sorted" and the other vectors are reshuffled
00521    // in the same way as necessary to have the first vector sorted. 
00522    // I.e. the correlation between the elements is kept.
00523    UInt_t nArrays=v.size();
00524    Double_t temp;
00525    if (nArrays > 0) {
00526       UInt_t sizeofarray=v[0].size();
00527       for (UInt_t i=0; i<sizeofarray; i++) {
00528          for (UInt_t j=sizeofarray-1; j>i; j--) {
00529             if (v[0][j-1] < v[0][j]) {
00530                for (UInt_t k=0; k< nArrays; k++) {
00531                   temp = v[k][j-1]; v[k][j-1] = v[k][j]; v[k][j] = temp;
00532                }
00533                if (NULL != vs) {
00534                   TString temps = (*vs)[j-1]; (*vs)[j-1] = (*vs)[j]; (*vs)[j] = temps;
00535                }
00536             }
00537          }
00538       }
00539    }
00540 }
00541 
00542 //_______________________________________________________________________
00543 Double_t TMVA::Tools::GetMutualInformation( const TH2F& h_ )
00544 {
00545    // Mutual Information method for non-linear correlations estimates in 2D histogram
00546    // Author: Moritz Backes, Geneva (2009)
00547 
00548    Double_t hi = h_.Integral();
00549    if (hi == 0) return -1; 
00550 
00551    // copy histogram and rebin to speed up procedure
00552    TH2F h( h_ );
00553    h.RebinX(2);
00554    h.RebinY(2);
00555    
00556    Double_t mutualInfo = 0.;
00557    Int_t maxBinX = h.GetNbinsX();
00558    Int_t maxBinY = h.GetNbinsY();
00559    for (Int_t x = 1; x <= maxBinX; x++) {
00560       for (Int_t y = 1; y <= maxBinY; y++) {
00561          Double_t p_xy = h.GetBinContent(x,y)/hi;
00562          Double_t p_x  = h.Integral(x,x,1,maxBinY)/hi;
00563          Double_t p_y  = h.Integral(1,maxBinX,y,y)/hi;
00564          if (p_x > 0. && p_y > 0. && p_xy > 0.){
00565             mutualInfo += p_xy*TMath::Log(p_xy / (p_x * p_y));
00566          }
00567       }
00568    }
00569 
00570    return mutualInfo;
00571 }
00572 
00573 //_______________________________________________________________________
00574 Double_t TMVA::Tools::GetCorrelationRatio( const TH2F& h_ )
00575 {
00576    // Compute Correlation Ratio of 2D histogram to estimate functional dependency between two variables
00577    // Author: Moritz Backes, Geneva (2009)
00578 
00579    Double_t hi = h_.Integral();
00580    if (hi == 0.) return -1; 
00581 
00582    // copy histogram and rebin to speed up procedure
00583    TH2F h( h_ );
00584    h.RebinX(2);
00585    h.RebinY(2);
00586 
00587    Double_t corrRatio = 0.;    
00588    Double_t y_mean = h.ProjectionY()->GetMean();
00589    for (Int_t ix=1; ix<=h.GetNbinsX(); ix++) {
00590       corrRatio += (h.Integral(ix,ix,1,h.GetNbinsY())/hi)*pow((GetYMean_binX(h,ix)-y_mean),2);
00591    }
00592    corrRatio /= pow(h.ProjectionY()->GetRMS(),2);
00593    return corrRatio;
00594 }
00595 
00596 //_______________________________________________________________________
00597 Double_t TMVA::Tools::GetYMean_binX( const TH2& h, Int_t bin_x )
00598 {
00599    // Compute the mean in Y for a given bin X of a 2D histogram
00600  
00601    if (h.Integral(bin_x,bin_x,1,h.GetNbinsY()) == 0.) {return 0;}
00602    Double_t y_bin_mean = 0.;
00603    TH1* py = h.ProjectionY();
00604    for (Int_t y = 1; y <= h.GetNbinsY(); y++){
00605       y_bin_mean += h.GetBinContent(bin_x,y)*py->GetBinCenter(y);
00606    }
00607    y_bin_mean /= h.Integral(bin_x,bin_x,1,h.GetNbinsY());
00608    return y_bin_mean;
00609 }
00610 
00611 //_______________________________________________________________________
00612 TH2F* TMVA::Tools::TransposeHist( const TH2F& h )
00613 {
00614    // Transpose quadratic histogram
00615 
00616    // sanity check
00617    if (h.GetNbinsX() != h.GetNbinsY()) {
00618       Log() << kFATAL << "<TransposeHist> cannot transpose non-quadratic histogram" << endl;
00619    }
00620    
00621    TH2F *transposedHisto = new TH2F( h ); 
00622    for (Int_t ix=1; ix <= h.GetNbinsX(); ix++){
00623       for (Int_t iy=1; iy <= h.GetNbinsY(); iy++){
00624          transposedHisto->SetBinContent(iy,ix,h.GetBinContent(ix,iy));
00625       }
00626    }
00627    return transposedHisto; // ownership returned
00628 }
00629 
00630 //_______________________________________________________________________
00631 Bool_t TMVA::Tools::CheckForSilentOption( const TString& cs ) const
00632 {
00633    // check for "silence" option in configuration option string
00634    Bool_t isSilent = kFALSE;
00635 
00636    TString s( cs ); 
00637    s.ToLower(); 
00638    s.ReplaceAll(" ","");
00639    if (s.Contains("silent") && !s.Contains("silent=f")) {
00640       if (!s.Contains("!silent") || s.Contains("silent=t")) isSilent = kTRUE;
00641    }
00642 
00643    return isSilent;
00644 }
00645 
00646 //_______________________________________________________________________
00647 Bool_t TMVA::Tools::CheckForVerboseOption( const TString& cs ) const
00648 {
00649    // check if verbosity "V" set in option
00650    Bool_t isVerbose = kFALSE;
00651 
00652    TString s( cs ); 
00653    s.ToLower(); 
00654    s.ReplaceAll(" ","");
00655    std::vector<TString> v = SplitString( s, ':' );
00656    for (std::vector<TString>::iterator it = v.begin(); it != v.end(); it++) {
00657       if ((*it == "v" || *it == "verbose") && !it->Contains("!")) isVerbose = kTRUE;
00658    }
00659 
00660    return isVerbose;
00661 }
00662 
00663 //_______________________________________________________________________
00664 void TMVA::Tools::UsefulSortDescending( std::vector<Double_t>& v )
00665 {
00666    // sort vector
00667    vector< vector<Double_t> > vtemp;
00668    vtemp.push_back(v);
00669    UsefulSortDescending(vtemp);
00670    v = vtemp[0];
00671 }
00672 
00673 //_______________________________________________________________________
00674 void TMVA::Tools::UsefulSortAscending( std::vector<Double_t>& v )
00675 {
00676    // sort vector
00677    vector<vector<Double_t> > vtemp;
00678    vtemp.push_back(v);
00679    UsefulSortAscending(vtemp);
00680    v = vtemp[0];
00681 }
00682 
00683 //_______________________________________________________________________
00684 Int_t TMVA::Tools::GetIndexMaxElement( std::vector<Double_t>& v )
00685 {
00686    // find index of maximum entry in vector
00687    if (v.size()==0) return -1;
00688 
00689    Int_t pos=0; Double_t mx=v[0];
00690    for (UInt_t i=0; i<v.size(); i++){
00691       if (v[i] > mx){
00692          mx=v[i];
00693          pos=i;
00694       }
00695    }
00696    return pos;
00697 }
00698 
00699 //_______________________________________________________________________
00700 Int_t TMVA::Tools::GetIndexMinElement( std::vector<Double_t>& v )
00701 {
00702    // find index of minimum entry in vector
00703    if (v.size()==0) return -1;
00704 
00705    Int_t pos=0; Double_t mn=v[0];
00706    for (UInt_t i=0; i<v.size(); i++){
00707       if (v[i] < mn){
00708          mn=v[i];
00709          pos=i;
00710       }
00711    }
00712    return pos;
00713 }
00714 
00715 
00716 //_______________________________________________________________________
00717 Bool_t TMVA::Tools::ContainsRegularExpression( const TString& s )  
00718 {
00719    // check if regular expression
00720    // helper function to search for "$!%^&()'<>?= " in a string
00721 
00722    Bool_t  regular = kFALSE;
00723    for (Int_t i = 0; i < Tools::fRegexp.Length(); i++) 
00724       if (s.Contains( Tools::fRegexp[i] )) { regular = kTRUE; break; }
00725 
00726    return regular;
00727 }
00728 
00729 //_______________________________________________________________________
00730 TString TMVA::Tools::ReplaceRegularExpressions( const TString& s, const TString& r )  
00731 {
00732    // replace regular expressions
00733    // helper function to remove all occurences "$!%^&()'<>?= " from a string
00734    // and replace all ::,$,*,/,+,- with _M_,_S_,_T_,_D_,_P_,_M_ respectively
00735 
00736    TString snew = s;
00737    for (Int_t i = 0; i < Tools::fRegexp.Length(); i++) 
00738       snew.ReplaceAll( Tools::fRegexp[i], r );
00739 
00740    snew.ReplaceAll( "::", r );
00741    snew.ReplaceAll( "$", "_S_" );
00742    snew.ReplaceAll( "&", "_A_" );
00743    snew.ReplaceAll( "%", "_MOD_" );
00744    snew.ReplaceAll( "|", "_O_" );
00745    snew.ReplaceAll( "*", "_T_" );
00746    snew.ReplaceAll( "/", "_D_" );
00747    snew.ReplaceAll( "+", "_P_" );
00748    snew.ReplaceAll( "-", "_M_" );
00749    snew.ReplaceAll( " ", "_" );
00750    snew.ReplaceAll( "[", "_" );
00751    snew.ReplaceAll( "]", "_" );
00752    snew.ReplaceAll( "=", "_E_" );
00753    snew.ReplaceAll( ">", "_GT_" );
00754    snew.ReplaceAll( "<", "_LT_" );
00755    snew.ReplaceAll( "(", "_" );
00756    snew.ReplaceAll( ")", "_" );
00757 
00758    return snew;
00759 }
00760 
00761 //_______________________________________________________________________
00762 const TString& TMVA::Tools::Color( const TString& c ) 
00763 {
00764    // human readable color strings
00765    static TString gClr_none         = "" ;
00766    static TString gClr_white        = "\033[1;37m";  // white
00767    static TString gClr_black        = "\033[30m";    // black
00768    static TString gClr_blue         = "\033[34m";    // blue
00769    static TString gClr_red          = "\033[1;31m" ; // red
00770    static TString gClr_yellow       = "\033[1;33m";  // yellow
00771    static TString gClr_darkred      = "\033[31m";    // dark red
00772    static TString gClr_darkgreen    = "\033[32m";    // dark green
00773    static TString gClr_darkyellow   = "\033[33m";    // dark yellow
00774                                     
00775    static TString gClr_bold         = "\033[1m"    ; // bold 
00776    static TString gClr_black_b      = "\033[30m"   ; // bold black
00777    static TString gClr_lblue_b      = "\033[1;34m" ; // bold light blue
00778    static TString gClr_cyan_b       = "\033[0;36m" ; // bold cyan
00779    static TString gClr_lgreen_b     = "\033[1;32m";  // bold light green
00780                                     
00781    static TString gClr_blue_bg      = "\033[44m";    // blue background
00782    static TString gClr_red_bg       = "\033[1;41m";  // white on red background
00783    static TString gClr_whiteonblue  = "\033[1;44m";  // white on blue background
00784    static TString gClr_whiteongreen = "\033[1;42m";  // white on green background
00785    static TString gClr_grey_bg      = "\033[47m";    // grey background
00786 
00787    static TString gClr_reset  = "\033[0m";     // reset
00788 
00789    if (!gConfig().UseColor()) return gClr_none;
00790 
00791    if (c == "white" )         return gClr_white; 
00792    if (c == "blue"  )         return gClr_blue; 
00793    if (c == "black"  )        return gClr_black; 
00794    if (c == "lightblue")      return gClr_cyan_b;
00795    if (c == "yellow")         return gClr_yellow; 
00796    if (c == "red"   )         return gClr_red; 
00797    if (c == "dred"  )         return gClr_darkred; 
00798    if (c == "dgreen")         return gClr_darkgreen; 
00799    if (c == "lgreenb")        return gClr_lgreen_b;
00800    if (c == "dyellow")        return gClr_darkyellow; 
00801 
00802    if (c == "bold")           return gClr_bold; 
00803    if (c == "bblack")         return gClr_black_b; 
00804 
00805    if (c == "blue_bgd")       return gClr_blue_bg; 
00806    if (c == "red_bgd" )       return gClr_red_bg; 
00807  
00808    if (c == "white_on_blue" ) return gClr_whiteonblue; 
00809    if (c == "white_on_green") return gClr_whiteongreen; 
00810 
00811    if (c == "reset") return gClr_reset; 
00812 
00813    cout << "Unknown color " << c << endl;
00814    exit(1);
00815 
00816    return gClr_none;
00817 }
00818 
00819 //_______________________________________________________________________
00820 void TMVA::Tools::FormattedOutput( const std::vector<Double_t>& values, const std::vector<TString>& V, 
00821                                    const TString titleVars, const TString titleValues, MsgLogger& logger,
00822                                    TString format )
00823 {
00824    // formatted output of simple table
00825 
00826    // sanity check
00827    UInt_t nvar = V.size();
00828    if ((UInt_t)values.size() != nvar) {
00829       logger << kFATAL << "<FormattedOutput> fatal error with dimensions: " 
00830              << values.size() << " OR " << " != " << nvar << Endl;
00831    }
00832 
00833    // find maximum length in V (and column title)
00834    UInt_t maxL = 7;
00835    std::vector<UInt_t> vLengths;
00836    for (UInt_t ivar=0; ivar<nvar; ivar++) maxL = TMath::Max( (UInt_t)V[ivar].Length(), maxL );
00837    maxL = TMath::Max( (UInt_t)titleVars.Length(), maxL );
00838 
00839    // column length
00840    UInt_t maxV = 7;
00841    maxV = TMath::Max( (UInt_t)titleValues.Length() + 1, maxL );
00842 
00843    // full column length
00844    UInt_t clen = maxL + maxV + 3;
00845 
00846    // bar line
00847    for (UInt_t i=0; i<clen; i++) logger << "-";
00848    logger << Endl;
00849 
00850    // title bar   
00851    logger << setw(maxL) << titleVars << ":";
00852    logger << setw(maxV+1) << titleValues << ":";
00853    logger << Endl;
00854    for (UInt_t i=0; i<clen; i++) logger << "-";
00855    logger << Endl;
00856 
00857    // the numbers
00858    for (UInt_t irow=0; irow<nvar; irow++) {
00859       logger << setw(maxL) << V[irow] << ":";
00860       logger << setw(maxV+1) << Form( format.Data(), values[irow] );
00861       logger << Endl;
00862    }
00863 
00864    // bar line
00865    for (UInt_t i=0; i<clen; i++) logger << "-";
00866    logger << Endl;
00867 }
00868 
00869 //_______________________________________________________________________
00870 void TMVA::Tools::FormattedOutput( const TMatrixD& M, const std::vector<TString>& V, MsgLogger& logger )
00871 {
00872    // formatted output of matrix (with labels)
00873 
00874    // sanity check: matrix must be quadratic
00875    UInt_t nvar = V.size();
00876    if ((UInt_t)M.GetNcols() != nvar || (UInt_t)M.GetNrows() != nvar) {
00877       logger << kFATAL << "<FormattedOutput> fatal error with dimensions: " 
00878              << M.GetNcols() << " OR " << M.GetNrows() << " != " << nvar << " ==> abort" << Endl;
00879    }
00880 
00881    // get length of each variable, and maximum length  
00882    UInt_t minL = 7;
00883    UInt_t maxL = minL;
00884    std::vector<UInt_t> vLengths;
00885    for (UInt_t ivar=0; ivar<nvar; ivar++) {
00886       vLengths.push_back(TMath::Max( (UInt_t)V[ivar].Length(), minL ));
00887       maxL = TMath::Max( vLengths.back(), maxL );
00888    }
00889    
00890    // count column length
00891    UInt_t clen = maxL+1;
00892    for (UInt_t icol=0; icol<nvar; icol++) clen += vLengths[icol]+1;
00893 
00894    // bar line
00895    for (UInt_t i=0; i<clen; i++) logger << "-";
00896    logger << Endl;
00897 
00898    // title bar   
00899    logger << setw(maxL+1) << " ";
00900    for (UInt_t icol=0; icol<nvar; icol++) logger << setw(vLengths[icol]+1) << V[icol];
00901    logger << Endl;
00902 
00903    // the numbers
00904    for (UInt_t irow=0; irow<nvar; irow++) {
00905       logger << setw(maxL) << V[irow] << ":";
00906       for (UInt_t icol=0; icol<nvar; icol++) {
00907          logger << setw(vLengths[icol]+1) << Form( "%+1.3f", M(irow,icol) );
00908       }      
00909       logger << Endl;
00910    }
00911 
00912    // bar line
00913    for (UInt_t i=0; i<clen; i++) logger << "-";
00914    logger << Endl;
00915 }
00916 
00917 //_______________________________________________________________________
00918 void TMVA::Tools::FormattedOutput( const TMatrixD& M, 
00919                                    const std::vector<TString>& vert, const std::vector<TString>& horiz, 
00920                                    MsgLogger& logger )
00921 {
00922    // formatted output of matrix (with labels)
00923 
00924    // sanity check: matrix must be quadratic
00925    UInt_t nvvar = vert.size();   
00926    UInt_t nhvar = horiz.size();
00927 
00928    // get length of each variable, and maximum length  
00929    UInt_t minL = 7;
00930    UInt_t maxL = minL;
00931    std::vector<UInt_t> vLengths;
00932    for (UInt_t ivar=0; ivar<nvvar; ivar++) {
00933       vLengths.push_back(TMath::Max( (UInt_t)vert[ivar].Length(), minL ));
00934       maxL = TMath::Max( vLengths.back(), maxL );
00935    }
00936    
00937    // count column length
00938    UInt_t minLh = 7;
00939    UInt_t maxLh = minLh;
00940    std::vector<UInt_t> hLengths;
00941    for (UInt_t ivar=0; ivar<nhvar; ivar++) {
00942       hLengths.push_back(TMath::Max( (UInt_t)horiz[ivar].Length(), minL ));
00943       maxLh = TMath::Max( hLengths.back(), maxLh );
00944    }
00945 
00946    UInt_t clen = maxLh+1;
00947    for (UInt_t icol=0; icol<nhvar; icol++) clen += hLengths[icol]+1;
00948 
00949    // bar line
00950    for (UInt_t i=0; i<clen; i++) logger << "-";
00951    logger << Endl;
00952 
00953    // title bar   
00954    logger << setw(maxL+1) << " ";
00955    for (UInt_t icol=0; icol<nhvar; icol++) logger << setw(hLengths[icol]+1) << horiz[icol];
00956    logger << Endl;
00957 
00958    // the numbers
00959    for (UInt_t irow=0; irow<nvvar; irow++) {
00960       logger << setw(maxL) << vert[irow] << ":";
00961       for (UInt_t icol=0; icol<nhvar; icol++) {
00962          logger << setw(hLengths[icol]+1) << Form( "%+1.3f", M(irow,icol) );
00963       }      
00964       logger << Endl;
00965    }
00966 
00967    // bar line
00968    for (UInt_t i=0; i<clen; i++) logger << "-";
00969    logger << Endl;
00970 }
00971 
00972 //_______________________________________________________________________
00973 TString TMVA::Tools::GetXTitleWithUnit( const TString& title, const TString& unit )
00974 {
00975    // histogramming utility
00976    return ( unit == "" ? title : ( title + "  [" + unit + "]" ) );
00977 }
00978 
00979 //_______________________________________________________________________
00980 TString TMVA::Tools::GetYTitleWithUnit( const TH1& h, const TString& unit, Bool_t normalised )
00981 {
00982    // histogramming utility
00983    TString retval = ( normalised ? "(1/N) " : "" );
00984    retval += Form( "dN_{ }/^{ }%.3g %s", h.GetXaxis()->GetBinWidth(1), unit.Data() );
00985    return retval;
00986 }
00987 
00988 //_______________________________________________________________________
00989 void TMVA::Tools::WriteFloatArbitraryPrecision( Float_t val, ostream& os )
00990 {
00991    // writes a float value with the available precision to a stream
00992    os << val << " :: ";
00993    void * c = &val;
00994    for (int i=0; i<4; i++) {
00995       Int_t ic = *((char*)c+i)-'\0';
00996       if (ic<0) ic+=256;
00997       os << ic << " ";
00998    }
00999    os << ":: ";
01000 }
01001 
01002 //_______________________________________________________________________
01003 void TMVA::Tools::ReadFloatArbitraryPrecision( Float_t& val, istream& is )
01004 {
01005    // reads a float value with the available precision from a stream
01006    Float_t a = 0;
01007    is >> a;
01008    TString dn;
01009    is >> dn;
01010    Int_t c[4];
01011    void * ap = &a;
01012    for (int i=0; i<4; i++) {
01013       is >> c[i];
01014       *((char*)ap+i) = '\0'+c[i];
01015    }
01016    is >> dn;
01017    val = a;
01018 }
01019 
01020 
01021 // XML file reading/writing helper functions
01022 
01023 //_______________________________________________________________________
01024 Bool_t TMVA::Tools::HasAttr( void* node, const char* attrname )
01025 {
01026    // add attribute from xml
01027    return xmlengine().HasAttr(node, attrname);
01028 }
01029 
01030 //_______________________________________________________________________
01031 void TMVA::Tools::ReadAttr( void* node, const char* attrname, TString& value )
01032 {
01033    // add attribute from xml
01034    if(!HasAttr(node, attrname)) {
01035       const char * nodename = xmlengine().GetNodeName(node);
01036       Log() << kFATAL << "Trying to read non-existing attribute '" << attrname << "' from xml node '" << nodename << "'" << Endl;
01037    }
01038    const char* val = xmlengine().GetAttr(node, attrname);
01039    value = TString(val);
01040 }
01041 
01042 //_______________________________________________________________________
01043 void TMVA::Tools::AddAttr( void* node, const char* attrname, const char* value )
01044 {
01045    // add attribute to node
01046    if( node == 0 ) return;
01047    gTools().xmlengine().NewAttr(node, 0, attrname, value );
01048 }
01049 
01050 //_______________________________________________________________________
01051 void* TMVA::Tools::AddChild( void* parent, const char* childname, const char* content, bool isRootNode ) {
01052    if( !isRootNode && parent == 0 ) return 0;
01053    return gTools().xmlengine().NewChild(parent, 0, childname, content);
01054 }
01055 
01056 //_______________________________________________________________________
01057 Bool_t TMVA::Tools::AddComment( void* node, const char* comment ) {
01058    if( node == 0 ) return kFALSE;
01059    return gTools().xmlengine().AddComment(node, comment);
01060 }
01061  //_______________________________________________________________________
01062 void* TMVA::Tools::GetParent( void* child)
01063 {
01064    void* par = xmlengine().GetParent(child);
01065    
01066    return par;
01067 }
01068 //_______________________________________________________________________
01069 void* TMVA::Tools::GetChild( void* parent, const char* childname )
01070 {
01071    void* ch = xmlengine().GetChild(parent);
01072    if (childname != 0) {
01073       while (ch!=0 && strcmp(xmlengine().GetNodeName(ch),childname) != 0) ch = xmlengine().GetNext(ch);
01074    }
01075    return ch;
01076 }
01077 
01078 //_______________________________________________________________________
01079 void* TMVA::Tools::GetNextChild( void* prevchild, const char* childname )
01080 {
01081    // XML helpers
01082    void* ch = xmlengine().GetNext(prevchild);
01083    if (childname != 0) {
01084       while (ch!=0 && strcmp(xmlengine().GetNodeName(ch),childname)!=0) ch = xmlengine().GetNext(ch);
01085    }
01086    return ch;
01087 }
01088 
01089 //_______________________________________________________________________
01090 const char* TMVA::Tools::GetContent( void* node )
01091 {
01092    // XML helpers
01093    return xmlengine().GetNodeContent(node);
01094 }
01095 
01096 //_______________________________________________________________________
01097 const char* TMVA::Tools::GetName( void* node )
01098 {
01099    // XML helpers
01100    return xmlengine().GetNodeName(node);
01101 }
01102 
01103 //_______________________________________________________________________
01104 Bool_t TMVA::Tools::AddRawLine( void* node, const char * raw )
01105 {
01106    // XML helpers
01107    return xmlengine().AddRawLine( node, raw );
01108 }
01109 
01110 //_______________________________________________________________________
01111 std::vector<TString> TMVA::Tools::SplitString(const TString& theOpt, const char separator ) const
01112 {
01113    // splits the option string at 'separator' and fills the list
01114    // 'splitV' with the primitive strings
01115    std::vector<TString> splitV;
01116    TString splitOpt(theOpt);
01117    splitOpt.ReplaceAll("\n"," ");
01118    splitOpt = splitOpt.Strip(TString::kBoth,separator);
01119    while (splitOpt.Length()>0) {
01120       if ( !splitOpt.Contains(separator) ) {
01121          splitV.push_back(splitOpt);
01122          break;
01123       }
01124       else {
01125          TString toSave = splitOpt(0,splitOpt.First(separator));
01126          splitV.push_back(toSave);
01127          splitOpt = splitOpt(splitOpt.First(separator),splitOpt.Length());
01128       }
01129       splitOpt = splitOpt.Strip(TString::kLeading,separator);
01130    }
01131    return splitV;
01132 }
01133 
01134 //_______________________________________________________________________
01135 TString TMVA::Tools::StringFromInt( Long_t i ) 
01136 {
01137    // string tools
01138    std::stringstream s;
01139    s << i;
01140    return TString(s.str().c_str());
01141 }
01142 
01143 //_______________________________________________________________________
01144 TString TMVA::Tools::StringFromDouble( Double_t d ) 
01145 {
01146    // string tools
01147    std::stringstream s;
01148    s << d;
01149    return TString(s.str().c_str());
01150 }
01151 
01152 //_______________________________________________________________________
01153 void TMVA::Tools::WriteTMatrixDToXML(void* node, const char* name, TMatrixD* mat)
01154 {
01155    // XML helpers
01156    void* matnode = xmlengine().NewChild(node, 0, name);
01157    xmlengine().NewAttr(matnode,0,"Rows", StringFromInt(mat->GetNrows()) );
01158    xmlengine().NewAttr(matnode,0,"Columns", StringFromInt(mat->GetNcols()) );
01159    std::stringstream s;
01160    for (Int_t row = 0; row<mat->GetNrows(); row++) {
01161       for (Int_t col = 0; col<mat->GetNcols(); col++) {
01162          s << (*mat)[row][col] << " ";
01163       }
01164    }
01165    xmlengine().AddRawLine( matnode, s.str().c_str() );
01166 }
01167 
01168 //_______________________________________________________________________
01169 void TMVA::Tools::WriteTVectorDToXML(void* node, const char* name, TVectorD* vec)
01170 {
01171    TMatrixD mat(1,vec->GetNoElements(),&((*vec)[0]));
01172    WriteTMatrixDToXML(node, name, &mat);
01173 }
01174 
01175 //_______________________________________________________________________
01176 void TMVA::Tools::ReadTVectorDFromXML(void* node, const char* name, TVectorD* vec)
01177 {
01178    TMatrixD mat(1,vec->GetNoElements(),&((*vec)[0]));
01179    ReadTMatrixDFromXML(node,name,&mat);
01180    for (int i=0;i<vec->GetNoElements();++i){
01181       (*vec)[i]=mat[0][i];
01182    }
01183 }
01184 
01185 //_______________________________________________________________________
01186 void TMVA::Tools::ReadTMatrixDFromXML(void* node, const char* name, TMatrixD* mat)
01187 {
01188    if (strcmp(xmlengine().GetNodeName(node),name)!=0){
01189       Log() << kWARNING << "Possible Error: Name of matrix in weight file"
01190             << " does not match name of matrix passed as argument!" << Endl;
01191    }
01192    Int_t nrows, ncols;
01193    ReadAttr(node, "Rows", nrows);
01194    ReadAttr(node, "Columns", ncols);
01195    if (mat->GetNrows() != nrows || mat->GetNcols() != ncols){
01196       Log() << kWARNING << "Possible Error: Dimension of matrix in weight file"
01197             << " does not match dimension of matrix passed as argument!" << Endl;
01198       mat->ResizeTo(nrows,ncols);
01199    }
01200    const char* content = xmlengine().GetNodeContent(node);
01201    std::stringstream s(content);
01202    for (Int_t row = 0; row<nrows; row++) {
01203       for (Int_t col = 0; col<ncols; col++) {
01204          s >> (*mat)[row][col];
01205       }
01206    }
01207 }
01208 
01209 //_______________________________________________________________________
01210 void TMVA::Tools::TMVAWelcomeMessage()
01211 {
01212    // direct output, eg, when starting ROOT session -> no use of Logger here
01213    cout << endl;
01214    cout << Color("bold") << "TMVA -- Toolkit for Multivariate Data Analysis" << Color("reset") << endl;
01215    cout << "        " << "Version " << TMVA_RELEASE << ", " << TMVA_RELEASE_DATE << endl;
01216    cout << "        " << "Copyright (C) 2005-2010 CERN, MPI-K Heidelberg, Us of Bonn and Victoria" << endl;
01217    cout << "        " << "Home page:     http://tmva.sf.net" << endl;
01218    cout << "        " << "Citation info: http://tmva.sf.net/citeTMVA.html" << endl;
01219    cout << "        " << "License:       http://tmva.sf.net/LICENSE" << endl << endl;
01220 }
01221 
01222 //_______________________________________________________________________
01223 void TMVA::Tools::TMVAVersionMessage( MsgLogger& logger )
01224 {
01225    // prints the TMVA release number and date
01226    logger << "___________TMVA Version " << TMVA_RELEASE << ", " << TMVA_RELEASE_DATE 
01227           << "" << Endl;
01228 }
01229 
01230 //_______________________________________________________________________
01231 void TMVA::Tools::ROOTVersionMessage( MsgLogger& logger )
01232 {
01233    // prints the ROOT release number and date
01234    static const char *months[] = { "Jan","Feb","Mar","Apr","May",
01235                                    "Jun","Jul","Aug","Sep","Oct",
01236                                    "Nov","Dec" };
01237    Int_t   idatqq = gROOT->GetVersionDate();   
01238    Int_t   iday   = idatqq%100;
01239    Int_t   imonth = (idatqq/100)%100;
01240    Int_t   iyear  = (idatqq/10000);
01241    TString versionDate = Form("%s %d, %4d",months[imonth-1],iday,iyear);
01242 
01243    logger << "You are running ROOT Version: " << gROOT->GetVersion() << ", " << versionDate << Endl;
01244 }
01245 
01246 //_______________________________________________________________________
01247 void TMVA::Tools::TMVAWelcomeMessage( MsgLogger& logger, EWelcomeMessage msgType )
01248 {
01249    // various kinds of welcome messages
01250    // ASCII text generated by this site: http://www.network-science.de/ascii/
01251 
01252    switch (msgType) {
01253 
01254    case kStandardWelcomeMsg:
01255       logger << Color("white") << "TMVA -- Toolkit for Multivariate Analysis" << Color("reset") << Endl;
01256       logger << "Copyright (C) 2005-2006 CERN, LAPP & MPI-K Heidelberg and Victoria U." << Endl;
01257       logger << "Home page http://tmva.sourceforge.net" << Endl;
01258       logger << "All rights reserved, please read http://tmva.sf.net/license.txt" << Endl << Endl;
01259       break;
01260 
01261    case kIsometricWelcomeMsg:
01262       logger << "   ___           ___           ___           ___      " << Endl;
01263       logger << "  /\\  \\         /\\__\\         /\\__\\         /\\  \\     " << Endl;
01264       logger << "  \\:\\  \\       /::|  |       /:/  /        /::\\  \\    " << Endl;
01265       logger << "   \\:\\  \\     /:|:|  |      /:/  /        /:/\\:\\  \\   " << Endl;
01266       logger << "   /::\\  \\   /:/|:|__|__   /:/__/  ___   /::\\~\\:\\  \\  " << Endl;
01267       logger << "  /:/\\:\\__\\ /:/ |::::\\__\\  |:|  | /\\__\\ /:/\\:\\ \\:\\__\\ " << Endl;
01268       logger << " /:/  \\/__/ \\/__/~~/:/  /  |:|  |/:/  / \\/__\\:\\/:/  / " << Endl;
01269       logger << "/:/  /            /:/  /   |:|__/:/  /       \\::/  /  " << Endl;
01270       logger << "\\/__/            /:/  /     \\::::/__/        /:/  /   " << Endl;
01271       logger << "                /:/  /       ~~~~           /:/  /    " << Endl;
01272       logger << "                \\/__/                       \\/__/     " << Endl << Endl;
01273       break;
01274 
01275    case kBlockWelcomeMsg:
01276       logger << Endl;
01277       logger << "_|_|_|_|_|  _|      _|  _|      _|    _|_|    " << Endl;
01278       logger << "    _|      _|_|  _|_|  _|      _|  _|    _|  " << Endl;
01279       logger << "    _|      _|  _|  _|  _|      _|  _|_|_|_|  " << Endl;
01280       logger << "    _|      _|      _|    _|  _|    _|    _|  " << Endl;
01281       logger << "    _|      _|      _|      _|      _|    _|  " << Endl << Endl;
01282       break;
01283 
01284    case kLeanWelcomeMsg:
01285       logger << Endl;
01286       logger << "_/_/_/_/_/  _/      _/  _/      _/    _/_/   " << Endl;
01287       logger << "   _/      _/_/  _/_/  _/      _/  _/    _/  " << Endl;
01288       logger << "  _/      _/  _/  _/  _/      _/  _/_/_/_/   " << Endl;
01289       logger << " _/      _/      _/    _/  _/    _/    _/    " << Endl;
01290       logger << "_/      _/      _/      _/      _/    _/     " << Endl << Endl;
01291       break;
01292 
01293    case kLogoWelcomeMsg:
01294       logger << Endl;
01295       logger << "_/_/_/_/_/ _|      _|  _|      _|    _|_|   " << Endl;
01296       logger << "   _/      _|_|  _|_|  _|      _|  _|    _| " << Endl;
01297       logger << "  _/       _|  _|  _|  _|      _|  _|_|_|_| " << Endl;
01298       logger << " _/        _|      _|    _|  _|    _|    _| " << Endl;
01299       logger << "_/         _|      _|      _|      _|    _| " << Endl << Endl;
01300       break;
01301 
01302    case kSmall1WelcomeMsg:
01303       logger << " _____ __  ____   ___   " << Endl;
01304       logger << "|_   _|  \\/  \\ \\ / /_\\  " << Endl;
01305       logger << "  | | | |\\/| |\\ V / _ \\ " << Endl;
01306       logger << "  |_| |_|  |_| \\_/_/ \\_\\" << Endl << Endl;
01307       break;
01308 
01309    case kSmall2WelcomeMsg:
01310       logger << " _____ __  ____     ___     " << Endl;
01311       logger << "|_   _|  \\/  \\ \\   / / \\    " << Endl;
01312       logger << "  | | | |\\/| |\\ \\ / / _ \\   " << Endl;
01313       logger << "  | | | |  | | \\ V / ___ \\  " << Endl;
01314       logger << "  |_| |_|  |_|  \\_/_/   \\_\\ " << Endl << Endl;
01315       break;
01316 
01317    case kOriginalWelcomeMsgColor:
01318       logger << kINFO << "" << Color("red") 
01319              << "_______________________________________" << Color("reset") << Endl;
01320       logger << kINFO << "" << Color("blue")
01321              << Color("red_bgd") << Color("bwhite") << " // " << Color("reset")
01322              << Color("white") << Color("blue_bgd") 
01323              << "|\\  /|| \\  //  /\\\\\\\\\\\\\\\\\\\\\\\\ \\ \\ \\ " << Color("reset") << Endl;
01324       logger << kINFO << ""<< Color("blue")
01325              << Color("red_bgd") << Color("white") << "//  " << Color("reset")
01326              << Color("white") << Color("blue_bgd") 
01327              << "| \\/ ||  \\//  /--\\\\\\\\\\\\\\\\\\\\\\\\ \\ \\ \\" << Color("reset") << Endl;
01328       break;
01329       
01330    case kOriginalWelcomeMsgBW:
01331       logger << kINFO << "" 
01332              << "_______________________________________" << Endl;
01333       logger << kINFO << " // "
01334              << "|\\  /|| \\  //  /\\\\\\\\\\\\\\\\\\\\\\\\ \\ \\ \\ " << Endl;
01335       logger << kINFO << "//  " 
01336              << "| \\/ ||  \\//  /--\\\\\\\\\\\\\\\\\\\\\\\\ \\ \\ \\" << Endl;
01337       break;
01338       
01339    default:
01340       logger << kFATAL << "unknown message type: " << msgType << Endl;
01341    }
01342 }
01343 
01344 //_______________________________________________________________________
01345 void TMVA::Tools::TMVACitation( MsgLogger& logger, ECitation citType )
01346 {
01347    // kinds of TMVA citation
01348 
01349    switch (citType) {
01350 
01351    case kPlainText:
01352       logger << "A. Hoecker, P. Speckmayer, J. Stelzer, J. Therhaag, E. von Toerne, H. Voss" << Endl;
01353       logger << "\"TMVA - Toolkit for Multivariate Data Analysis\" PoS ACAT:040,2007. e-Print: physics/0703039" << Endl;
01354       break;
01355 
01356    case kBibTeX:
01357       logger << "@Article{TMVA2007," << Endl;
01358       logger << "     author    = \"Hoecker, Andreas and Speckmayer, Peter and Stelzer, Joerg " << Endl;
01359       logger << "                   and Therhaag, Jan and von Toerne, Eckhard and Voss, Helge\"," << Endl;
01360       logger << "     title     = \"{TMVA: Toolkit for multivariate data analysis}\"," << Endl;
01361       logger << "     journal   = \"PoS\"," << Endl;
01362       logger << "     volume    = \"ACAT\"," << Endl;
01363       logger << "     year      = \"2007\"," << Endl;
01364       logger << "     pages     = \"040\"," << Endl;
01365       logger << "     eprint    = \"physics/0703039\"," << Endl;
01366       logger << "     archivePrefix = \"arXiv\"," << Endl;
01367       logger << "     SLACcitation  = \"%%CITATION = PHYSICS/0703039;%%\"" << Endl;
01368       logger << "}" << Endl;
01369       break;
01370 
01371    case kLaTeX:
01372       logger << "%\\cite{TMVA2007}" << Endl;
01373       logger << "\bibitem{TMVA2007}" << Endl;
01374       logger << "  A.~Hoecker, P.~Speckmayer, J.~Stelzer, J.~Therhaag, E.~von Toerne, H.~Voss" << Endl;
01375       logger << "  %``TMVA: Toolkit for multivariate data analysis,''" << Endl;
01376       logger << "  PoS A {\\bf CAT} (2007) 040" << Endl;
01377       logger << "  [arXiv:physics/0703039]." << Endl;
01378       logger << "  %%CITATION = POSCI,ACAT,040;%%" << Endl;
01379       break;
01380 
01381    case kHtmlLink:
01382       logger << kINFO << "  " << Endl;
01383       logger << kINFO << gTools().Color("bold") 
01384          << "Thank you for using TMVA!" << gTools().Color("reset") << Endl;
01385       logger << kINFO << gTools().Color("bold") 
01386              << "For citation information, please visit: http://tmva.sf.net/citeTMVA.html"
01387              << gTools().Color("reset") << Endl; 
01388    }
01389 }
01390 
01391 //_______________________________________________________________________
01392 Bool_t TMVA::Tools::HistoHasEquidistantBins(const TH1& h)
01393 {
01394    return !(h.GetXaxis()->GetXbins()->fN);
01395 }
01396 
01397 //_______________________________________________________________________
01398 std::vector<TMatrixDSym*>*
01399 TMVA::Tools::CalcCovarianceMatrices( const std::vector<Event*>& events, Int_t maxCls )
01400 {
01401    // compute covariance matrices
01402 
01403    if (events.size() == 0) return 0;
01404 
01405 
01406    UInt_t nvar = events.at(0)->GetNVariables(), ivar = 0, jvar = 0;
01407 
01408    // init matrices
01409    Int_t matNum = maxCls;
01410    if (maxCls > 1 ) matNum++; // if more than one classes, then produce one matrix for all events as well (beside the matrices for each class)
01411 
01412    std::vector<TVectorD*>* vec = new std::vector<TVectorD*>(matNum);
01413    std::vector<TMatrixD*>* mat2 = new std::vector<TMatrixD*>(matNum);
01414    std::vector<Double_t> count(matNum);
01415    count.assign(matNum,0);
01416 
01417    Int_t cls = 0;
01418    TVectorD* v;
01419    TMatrixD* m;
01420    for (cls = 0; cls < matNum ; cls++) {
01421       vec->at(cls) = new TVectorD(nvar);
01422       mat2->at(cls) = new TMatrixD(nvar,nvar);
01423       v = vec->at(cls);
01424       m = mat2->at(cls);
01425 
01426       for (ivar=0; ivar<nvar; ivar++) {
01427          (*v)(ivar) = 0;
01428          for (jvar=0; jvar<nvar; jvar++) {
01429             (*m)(ivar, jvar) = 0;
01430          }
01431       }
01432    }
01433 
01434    // perform event loop
01435    for (UInt_t i=0; i<events.size(); i++) {
01436 
01437       // fill the event
01438       Event * ev = events[i];
01439       cls = ev->GetClass();
01440       Double_t weight = ev->GetWeight();
01441        
01442       if (maxCls > 1) {
01443          v = vec->at(matNum-1);
01444          m = mat2->at(matNum-1);
01445 
01446          count.at(matNum-1)+=weight; // count used events
01447          for (ivar=0; ivar<nvar; ivar++) {
01448 
01449             Double_t xi = ev->GetValue(ivar);
01450             (*v)(ivar) += xi*weight;
01451             (*m)(ivar, ivar) += (xi*xi*weight);
01452 
01453             for (jvar=ivar+1; jvar<nvar; jvar++) {
01454                Double_t xj = ev->GetValue(jvar);
01455                (*m)(ivar, jvar) += (xi*xj*weight);
01456                (*m)(jvar, ivar) = (*m)(ivar, jvar); // symmetric matrix
01457             }
01458          }
01459       }
01460 
01461       count.at(cls)+=weight; // count used events
01462       v = vec->at(cls);
01463       m = mat2->at(cls);
01464       for (ivar=0; ivar<nvar; ivar++) {
01465          Double_t xi = ev->GetValue(ivar);
01466          (*v)(ivar) += xi*weight;
01467          (*m)(ivar, ivar) += (xi*xi*weight);
01468 
01469          for (jvar=ivar+1; jvar<nvar; jvar++) {
01470             Double_t xj = ev->GetValue(jvar);
01471             (*m)(ivar, jvar) += (xi*xj*weight);
01472             (*m)(jvar, ivar) = (*m)(ivar, jvar); // symmetric matrix
01473          }
01474       }
01475    }
01476 
01477    // variance-covariance
01478    std::vector<TMatrixDSym*>* mat = new std::vector<TMatrixDSym*>(matNum);
01479    for (cls = 0; cls < matNum; cls++) {
01480       v = vec->at(cls);
01481       m = mat2->at(cls);
01482 
01483       mat->at(cls) = new TMatrixDSym(nvar);
01484 
01485       Double_t n = count.at(cls);
01486       for (ivar=0; ivar<nvar; ivar++) {
01487          for (jvar=0; jvar<nvar; jvar++) {
01488             (*(mat->at(cls)))(ivar, jvar) = (*m)(ivar, jvar)/n - (*v)(ivar)*(*v)(jvar)/(n*n);
01489          }
01490       }
01491       delete v;
01492       delete m;
01493    }
01494 
01495    delete mat2;
01496    delete vec;
01497 
01498    return mat;
01499 }
01500 

Generated on Tue Jul 5 15:25:37 2011 for ROOT_528-00b_version by  doxygen 1.5.1