multidimSampling.C

Go to the documentation of this file.
00001 //+ example of sampling a multi-dim distribution using the DistSampler class
00002 //NOTE: This tutorial must be run with ACLIC
00003 //Author L. Moneta Dec 2010
00004 
00005 // function (a 4d gaussian)
00006 #include "TMath.h"
00007 #include "TF2.h"
00008 #include "TStopwatch.h"
00009 #include "Math/DistSampler.h"
00010 #include "Math/DistSamplerOptions.h"
00011 #include "Math/MinimizerOptions.h"
00012 #include "Math/Factory.h"
00013 
00014 #include "TKDTreeBinning.h"
00015 
00016 #include "TTree.h"
00017 #include "TFile.h"
00018 #include "TMatrixDSym.h"
00019 #include "TVectorD.h"
00020 #include "TCanvas.h"
00021 #include <cmath>
00022 
00023 // Gauss ND function
00024 // make a class in order to avoid constructing the 
00025 // matrices for every call
00026 // This however requires that  the code must be compiled with ACLIC
00027 
00028 bool debug = false;
00029 struct GausND { 
00030    
00031    TVectorD X; 
00032    TVectorD Mu; 
00033    TMatrixDSym CovMat; 
00034    
00035    GausND( int dim ) : 
00036       X(TVectorD(dim)), 
00037       Mu(TVectorD(dim)), 
00038       CovMat(TMatrixDSym(dim) )
00039    {}
00040    
00041    double operator() (double *x, double *p) { 
00042       // 4 parameters 
00043       int dim = X.GetNrows();
00044       int k = 0;
00045       for (int i = 0; i<dim; ++i) { X[i] = x[i] - p[k]; k++; } 
00046       for (int i = 0; i<dim; ++i) { 
00047          CovMat(i,i) = p[k]*p[k];
00048          k++;
00049       }
00050       for (int i = 0; i<dim; ++i) { 
00051          for (int j = i+1; j<dim; ++j) { 
00052             // p now are the correlations N(N-1)/2
00053                CovMat(i,j) = p[k]*sqrt(CovMat(i,i)*CovMat(j,j));            
00054                CovMat(j,i) = CovMat(i,j);
00055                k++;
00056          }
00057       }
00058       if (debug) { 
00059          X.Print();
00060          CovMat.Print();
00061       }
00062 
00063       double det = CovMat.Determinant(); 
00064       if (det <= 0) {
00065          Fatal("GausND","Determinant is <= 0 det = %f",det);
00066          CovMat.Print();
00067          return 0;  
00068       }
00069       double norm = std::pow( 2. * TMath::Pi(), dim/2) * sqrt(det);  
00070       // compute the gaussians
00071       CovMat.Invert();
00072       double fval  = std::exp( - 0.5 * CovMat.Similarity(X) )/ norm; 
00073 
00074       if (debug) { 
00075          std::cout << "det  " << det << std::endl; 
00076          std::cout << "norm " << norm << std::endl; 
00077          std::cout << "fval " << fval << std::endl; 
00078       }
00079 
00080       return fval;
00081    }
00082 };
00083 
00084 using namespace ROOT::Math; 
00085 
00086 void multidimSampling() { 
00087 
00088 #ifdef __CINT__
00089    std::cout << "DO NOT RUN WITH CINT:" << std::endl;
00090    std::cout << "we are using a custom function which requires" << std::endl;
00091    std::cout << "that this tutorial must be compiled with ACLIC" << std::endl;
00092   return;
00093 #endif
00094 
00095    const int N = 100000;
00096    //const int NBin = 1000;
00097    const int DIM = 4; 
00098 
00099    double xmin[] = {-10,-10,-10, -10};
00100    double xmax[] = { 10, 10, 10,  10};
00101    double par0[] = { 1., -1., 2, 0, // the gaussian mu
00102                      1, 2, 1, 3, // the sigma
00103                      0.5,0.,0.,0.,0.,0.8 };  // the correlation
00104    
00105    const int NPAR = DIM + DIM*(DIM+1)/2; // 14 in the 4d case
00106    // generate the sample 
00107    GausND gaus4d(4);
00108    TF1 * f = new TF1("functionND",gaus4d,0,1,14);
00109    f->SetParameters(par0);
00110 
00111    double x0[] = {0,0,0,0};
00112    // for debugging
00113    if (debug) f->EvalPar(x0,0);
00114    debug = false;
00115 
00116    TString name; 
00117    for (int i = 0; i < NPAR; ++i )  {
00118       if (i < DIM) f->SetParName(i, name.Format("mu_%d",i+1) );
00119       else if (i < 2*DIM) f->SetParName(i, name.Format("sig_%d",i-DIM+1) );
00120       else if (i < 2*DIM) f->SetParName(i, name.Format("sig_%d",i-2*DIM+1) );
00121    }
00122 
00123    //ROOT::Math::DistSamplerOptions::SetDefaultSampler("Foam");
00124    DistSampler * sampler = Factory::CreateDistSampler();
00125    if (sampler == 0) { 
00126       Info("multidimSampling","Default sampler %s is not available try with Foam ",
00127            ROOT::Math::DistSamplerOptions::DefaultSampler().c_str() );           
00128       ROOT::Math::DistSamplerOptions::SetDefaultSampler("Foam");
00129    }
00130    sampler = Factory::CreateDistSampler();
00131    if (sampler == 0) { 
00132       Error("multidimSampling","Foam sampler is not available - exit ");
00133       return;
00134    }
00135             
00136 
00137    sampler->SetFunction(*f,DIM);
00138    sampler->SetRange(xmin,xmax);
00139    bool ret = sampler->Init();
00140 
00141    std::vector<double> data1(DIM*N); 
00142    double v[DIM];
00143    TStopwatch w; 
00144 
00145    if (!ret) { 
00146       Error("Sampler::Init","Error initializing unuran sampler");
00147       return; 
00148    }
00149    // generate the data 
00150    w.Start(); 
00151    for (int i = 0; i < N; ++i) { 
00152       sampler->Sample(v);
00153       for (int j = 0; j < DIM; ++j) 
00154          data1[N*j + i]     = v[j];
00155    }
00156    w.Stop(); 
00157    w.Print(); 
00158 
00159    // fill tree with data 
00160    TFile * file = new TFile("multiDimSampling.root","RECREATE");
00161    double x[DIM];
00162    TTree * t1 = new TTree("t1","Tree from Unuran"); 
00163    t1->Branch("x",x,"x[4]/D");
00164    for (int i = 0; i < N; ++i) { 
00165       for (int j = 0; j < DIM; ++j) {
00166          x[j] = data1[i+N*j];
00167       }
00168       t1->Fill();
00169    }
00170 
00171 
00172    // try to fit the 2d data; 
00173    // GausND gaus2(2); 
00174    // TF2 * f2 = new TF2("f2",gaus2,-3,3,-3,3,5,"GausND");
00175    
00176    // f2->SetParameters(0,0,1,1,0);
00177    // f2->SetParLimits(4,-1,1);
00178    // t1->UnbinnedFit("f2","x[0]:x[1]");
00179 
00180    // plot the data
00181    
00182    t1->Draw("x[0]:x[1]:x[2]:x[3]","","candle");
00183    TCanvas * c2 = new TCanvas();
00184    c2->Divide(3,2);
00185    int ic=1;
00186    c2->cd(ic++);
00187    t1->Draw("x[0]:x[1]");
00188    c2->cd(ic++);
00189    t1->Draw("x[0]:x[2]");
00190    c2->cd(ic++);
00191    t1->Draw("x[0]:x[3]");
00192    c2->cd(ic++);
00193    t1->Draw("x[1]:x[2]");
00194    c2->cd(ic++);
00195    t1->Draw("x[1]:x[3]");
00196    c2->cd(ic++);
00197    t1->Draw("x[2]:x[3]");
00198    
00199 
00200    t1->Write();
00201    file->Close();
00202 
00203 
00204 }

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