]> git.uio.no Git - u/mrichter/AliRoot.git/blob - CORRFW/test/testCFContainers.C
Updates for real data analysis, trigger decisions
[u/mrichter/AliRoot.git] / CORRFW / test / testCFContainers.C
1 #include <Riostream.h>
2
3 extern TRandom *gRandom;
4 extern TBenchmark *gBenchmark;
5 extern TSystem *gSystem;
6
7 void testCFContainers(){
8
9   // simple example macros for usage of a N-dim container (AliCFContainer)
10   // handling a set of grids to accumulate data at different 
11   // selection steps, & derive efficiency (this is stored in AliCFEffGrid) 
12   // book, fill and draw some histos
13   // The efficiency is then used to correct the data (trivially self-correct, 
14   // in this example)
15
16   gROOT->SetStyle("Plain");
17   gStyle->SetPalette(1);
18   gStyle->SetOptStat(111110);
19   gStyle->SetPalette(1);
20   gStyle->SetCanvasColor(0);
21   gStyle->SetFrameFillColor(0);
22
23   gSystem->Load("libANALYSIS.so");
24   gSystem->Load("libCORRFW.so") ;
25  
26   //Setting up the container grid... 
27
28   const Int_t nstep=2; //number of selection steps  (just 2 in this ex)
29
30   const Int_t nvar=4; //number of variables on the grid:pt,vtx
31
32   const Int_t nbin1=6; //bins in pt
33   const Int_t nbin2=10; //bins in eta 
34   const Int_t nbin3=18; //bins in phi
35   const Int_t nbin4=10; //bins in vertex
36
37
38   //Flag the sel steps. In this example, we have two, may be any nstep
39   Int_t stepGen=0;
40   Int_t stepRec=1;
41
42   //the sensitive variables, their indeces
43   Int_t ipt =0;
44   Int_t ieta=1;
45   Int_t iphi=2;
46   Int_t ivtx=3;
47
48   //arrays for the number of bins in each dimension
49   const Int_t iBin[nvar] ={nbin1,nbin2,nbin3,nbin4};
50
51   //arrays for bin limits
52   Double_t binLim1[nbin1+1];
53   Double_t binLim2[nbin2+1];
54   Double_t binLim3[nbin3+1];
55   Double_t binLim4[nbin4+1];
56   
57   for(Int_t i=0;i<=nbin1;i++){
58     // pt [0,3] GeV/c
59     binLim1[i]=i*0.5; 
60   }
61
62   for(Int_t i=0;i<=nbin2;i++){
63     // eta [-1,1] 
64     binLim2[i]=i*0.2-1.;
65   }
66   for(Int_t i=0;i<=nbin3;i++){
67     //phi [0,360]
68     binLim3[i]=i*20.;
69   }
70   for(Int_t i=0;i<=nbin4;i++){
71     //vertex [-20,20] cm
72     binLim4[i]=i*4.-20.;
73   }
74   
75   
76   //the nstep grids "container" 
77   AliCFContainer *cont = new AliCFContainer("cont","example of  container",nstep,nvar,iBin);
78
79   //setting the bin limits
80   cont->SetBinLimits(ipt,binLim1);
81   cont->SetBinLimits(ieta,binLim2);
82   cont->SetBinLimits(iphi,binLim3);
83   cont->SetBinLimits(ivtx,binLim4);
84
85   //Start filling the mc and the data
86
87   //data sample (1M tracks)
88   Int_t nev=100000;
89   Int_t seed =1234;
90   gRandom->SetSeed(seed);
91   Double_t Value[nvar];
92   for(Int_t iev =0;iev<nev;iev++){
93     Float_t y=gRandom->Rndm();
94     Float_t pt=-TMath::Log(y)/0.5; //pt, exponential
95     Double_t eta=2.*gRandom->Rndm()-1.;//flat in eta
96     Double_t phi=360.*gRandom->Rndm(); //flat in phi
97     Float_t vtx=gRandom->Gaus( 0,5.);//gaussian in vertex
98     Value[ipt]=pt;
99     Value[ieta]=eta;
100     Value[iphi]=phi;
101     Value[ivtx]=vtx;    
102     cont->Fill(Value, stepGen); //fill the efficiency denominator, sel step=0
103     Float_t rndm=gRandom->Rndm();
104     //simulate 80% constant efficiency everywhere
105     if(rndm<0.8){
106       cont->Fill(Value,stepRec); //fill the efficiency denominator, sel step =1
107     }           
108   }   
109
110   //   Save it to a file
111    cont->Save("container.root");
112   //delete it
113    delete cont;
114
115 // Read the  container from file
116    TFile *file = new TFile("container.root");
117    AliCFContainer *data = (AliCFContainer*) (file->Get("cont"));
118
119   // Make some 1 & 2-D projections..
120   // pt and vertex, generator and reconstructed level
121 //   TCanvas *cmc =new TCanvas("cmc","The  distributions",0,300,900,900);
122 //   cmc->Divide(2,2);
123 //   cmc->cd(1);
124 //   TH1D *hpt1a = data->ShowProjection(ipt, stepGen);
125 //   hpt1a->SetMinimum(0.01);
126 //   hpt1a->Draw();
127 //   cmc->cd(2);
128 //   TH1D *hpt1b = data->ShowProjection(ipt, stepRec);
129 //   hpt1b->SetMinimum(0.01);
130 //   hpt1b->Draw();
131 //   cmc->cd(3);
132 //   TH2D *hptvtx1a = data->ShowProjection(ipt,ivtx, stepGen);
133 //   hptvtx1a->SetMinimum(0.01);
134 //   hptvtx1a->Draw("lego");
135 //   cmc->cd(4);
136 //   TH2D *hptvtx1b = data->ShowProjection(ipt,ivtx, stepRec);
137 //   hptvtx1b->SetMinimum(0.01);
138 //   hptvtx1b->Draw("lego");
139 //   cmc->Print("data.gif");
140
141  
142   //construct the efficiency grid from the data container 
143   AliCFEffGrid *eff = new AliCFEffGrid("eff"," The efficiency",*data);
144   eff->CalculateEfficiency(stepRec,stepGen); //eff= step1/step0
145
146   //The efficiency along pt and vertex, and 2-D projection
147   TCanvas *ceff =new TCanvas("ceff"," Efficiency",0,300,900,300);
148   ceff->Divide(3,1);
149   ceff->cd(1);
150   TH1D *hpt2a = eff->Project(ipt); //the efficiency vs pt
151   hpt2a->SetMinimum(0.01);
152   hpt2a->Draw();
153   ceff->cd(2);
154   TH1D *hvtx2a = eff->Project(ivtx); //the efficiency vs vtx
155   hvtx2a->SetMinimum(0.01);
156   hvtx2a->Draw();
157   ceff->cd(3);
158   TH2D *hptvtx2a = eff->Project(ipt,ivtx); //look at the numerator
159   hptvtx2a->SetMinimum(0.01);
160   hptvtx2a->SetMinimum(0.01);
161   hptvtx2a->Draw("lego");
162   
163   ceff->Print("eff.gif");
164
165   //get the corrected data grid  
166   AliCFDataGrid *corrdata = new AliCFDataGrid("corrdata","corrected data",*data);
167   //correct selection step "reconstructed"
168   corrdata->SetMeasured(stepRec); //set data to be corrected
169   corrdata->ApplyEffCorrection(*eff);//apply the correction for efficiency
170
171   //The observed data, the corrected ones and the "MC truth" distributions 
172   //vs pt and vtx
173   TCanvas *ccorrdata =new TCanvas("ccorrdata"," corrected data",0,300,900,900);
174   ccorrdata->Divide(2,2);
175   ccorrdata->cd(1);
176   TH1D *hpt3a = corrdata->GetData()->Project(ipt); //uncorrected data
177   hpt3a->SetMinimum(0.01);
178   hpt3a->Draw();
179   ccorrdata->cd(2);
180   TH1D *hpt3b = corrdata->Project(ipt); //corrected data
181   hpt3b->SetMinimum(0.01);
182   hpt3b->Draw();
183   ccorrdata->cd(3);
184   TH1D *hvtx3a = corrdata->GetData()->Project(ivtx); //uncorrected data
185   hvtx3a->SetMinimum(0.01);
186   hvtx3a->Draw();
187   ccorrdata->cd(4);
188   TH1D *hvtx3b = corrdata->Project(ivtx); //corrected data
189   hvtx3b->SetMinimum(0.01);
190   hvtx3b->Draw();
191   ccorrdata->Print("corrdata.gif");
192
193 }