rf401_importttreethx.C

Go to the documentation of this file.
00001 //////////////////////////////////////////////////////////////////////////
00002 //
00003 // 'DATA AND CATEGORIES' RooFit tutorial macro #401
00004 // 
00005 // Overview of advanced option for importing data from ROOT TTree and THx histograms
00006 // Basic import options are demonstrated in rf102_dataimport.C
00007 //
00008 //
00009 // 07/2008 - Wouter Verkerke 
00010 // 
00011 /////////////////////////////////////////////////////////////////////////
00012 
00013 #ifndef __CINT__
00014 #include "RooGlobalFunc.h"
00015 #endif
00016 #include "RooRealVar.h"
00017 #include "RooDataSet.h"
00018 #include "RooDataHist.h"
00019 #include "RooCategory.h"
00020 #include "RooGaussian.h"
00021 #include "RooConstVar.h"
00022 #include "TCanvas.h"
00023 #include "TAxis.h"
00024 #include "RooPlot.h"
00025 #include "TH1.h"
00026 #include "TTree.h"
00027 #include "TRandom.h"
00028 
00029 using namespace RooFit ;
00030 
00031 
00032 
00033 TH1* makeTH1(const char* name, Double_t mean, Double_t sigma) ;
00034 TTree* makeTTree() ;
00035 
00036 
00037 void rf401_importttreethx()
00038 {
00039   // I m p o r t  m u l t i p l e   T H 1   i n t o   a   R o o D a t a H i s t
00040   // --------------------------------------------------------------------------
00041 
00042   // Create thee ROOT TH1 histograms
00043   TH1* hh_1 = makeTH1("hh1",0,3) ;
00044   TH1* hh_2 = makeTH1("hh2",-3,1) ;
00045   TH1* hh_3 = makeTH1("hh3",+3,4) ;
00046 
00047   // Declare observable x
00048   RooRealVar x("x","x",-10,10) ;
00049 
00050   // Create category observable c that serves as index for the ROOT histograms
00051   RooCategory c("c","c") ;
00052   c.defineType("SampleA") ;
00053   c.defineType("SampleB") ;
00054   c.defineType("SampleC") ;
00055 
00056   // Create a binned dataset that imports contents of all TH1 mapped by index category c
00057   RooDataHist* dh = new RooDataHist("dh","dh",x,Index(c),Import("SampleA",*hh_1),Import("SampleB",*hh_2),Import("SampleC",*hh_3)) ;
00058   dh->Print() ;
00059 
00060   // Alternative constructor form for importing multiple histograms
00061   map<string,TH1*> hmap ;
00062   hmap["SampleA"] = hh_1 ;
00063   hmap["SampleB"] = hh_2 ;
00064   hmap["SampleC"] = hh_3 ;
00065   RooDataHist* dh2 = new RooDataHist("dh","dh",x,c,hmap) ;
00066   dh2->Print() ;
00067 
00068   
00069 
00070   // I m p o r t i n g   a   T T r e e   i n t o   a   R o o D a t a S e t   w i t h   c u t s 
00071   // -----------------------------------------------------------------------------------------
00072 
00073   TTree* tree = makeTTree() ;
00074 
00075   // Define observables y,z
00076   RooRealVar y("y","y",-10,10) ;
00077   RooRealVar z("z","z",-10,10) ;
00078 
00079   // Import only observables (y,z)
00080   RooDataSet ds("ds","ds",RooArgSet(x,y),Import(*tree)) ;
00081   ds.Print() ;
00082 
00083   // Import observables (x,y,z) but only event for which (y+z<0) is true
00084   RooDataSet ds2("ds2","ds2",RooArgSet(x,y,z),Import(*tree),Cut("y+z<0")) ;
00085   ds2.Print() ;
00086 
00087 
00088 
00089   // I m p o r t i n g   i n t e g e r   T T r e e   b r a n c h e s
00090   // ---------------------------------------------------------------
00091 
00092   // Import integer tree branch as RooRealVar
00093   RooRealVar i("i","i",0,5) ;
00094   RooDataSet ds3("ds3","ds3",RooArgSet(i,x),Import(*tree)) ;
00095   ds3.Print() ;
00096 
00097   // Define category i
00098   RooCategory icat("i","i") ;
00099   icat.defineType("State0",0) ;
00100   icat.defineType("State1",1) ;
00101 
00102   // Import integer tree branch as RooCategory (only events with i==0 and i==1
00103   // will be imported as those are the only defined states)
00104   RooDataSet ds4("ds4","ds4",RooArgSet(icat,x),Import(*tree)) ;
00105   ds4.Print() ;
00106 
00107 
00108 
00109   // I m p o r t  m u l t i p l e   R o o D a t a S e t s   i n t o   a   R o o D a t a S e t 
00110   // ----------------------------------------------------------------------------------------
00111 
00112   // Create three RooDataSets in (y,z)
00113   RooDataSet* dsA = (RooDataSet*) ds2.reduce(RooArgSet(x,y),"z<-5") ;
00114   RooDataSet* dsB = (RooDataSet*) ds2.reduce(RooArgSet(x,y),"abs(z)<5") ;
00115   RooDataSet* dsC = (RooDataSet*) ds2.reduce(RooArgSet(x,y),"z>5") ;
00116 
00117   // Create a dataset that imports contents of all the above datasets mapped by index category c
00118   RooDataSet* dsABC = new RooDataSet("dsABC","dsABC",RooArgSet(x,y),Index(c),Import("SampleA",*dsA),Import("SampleB",*dsB),Import("SampleC",*dsC)) ;
00119 
00120   dsABC->Print() ;
00121 
00122 }
00123 
00124 
00125 
00126 TH1* makeTH1(const char* name, Double_t mean, Double_t sigma) 
00127 {
00128   // Create ROOT TH1 filled with a Gaussian distribution
00129 
00130   TH1D* hh = new TH1D(name,name,100,-10,10) ;
00131   for (int i=0 ; i<1000 ; i++) {
00132     hh->Fill(gRandom->Gaus(mean,sigma)) ;
00133   }
00134   return hh ;
00135 }
00136 
00137 
00138 
00139 TTree* makeTTree() 
00140 {
00141   // Create ROOT TTree filled with a Gaussian distribution in x and a uniform distribution in y
00142 
00143   TTree* tree = new TTree("tree","tree") ;
00144   Double_t* px = new Double_t ;
00145   Double_t* py = new Double_t ;
00146   Double_t* pz = new Double_t ;
00147   Int_t*    pi = new Int_t ;
00148   tree->Branch("x",px,"x/D") ;
00149   tree->Branch("y",py,"y/D") ;
00150   tree->Branch("z",pz,"z/D") ;
00151   tree->Branch("i",pi,"i/I") ;
00152   for (int i=0 ; i<100 ; i++) {
00153     *px = gRandom->Gaus(0,3) ;
00154     *py = gRandom->Uniform()*30 - 15 ;
00155     *pz = gRandom->Gaus(0,5) ;
00156     *pi = i % 3 ;
00157     tree->Fill() ;
00158   }
00159   return tree ;
00160 }
00161 
00162 
00163 

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