]> git.uio.no Git - u/mrichter/AliRoot.git/blob - CORRFW/test/AliCFSingleTrackTaskCAF.C
QA histograms filling for missing steps
[u/mrichter/AliRoot.git] / CORRFW / test / AliCFSingleTrackTaskCAF.C
1 //DEFINITION OF A FEW CONSTANTS
2 const Double_t ymin  = -1.0 ;
3 const Double_t ymax  =  1.0 ;
4 const Double_t ptmin =  0.0 ;
5 const Double_t ptmax =  8.0 ;
6 const Int_t    mintrackrefsTPC = 2 ;
7 const Int_t    mintrackrefsITS = 3 ;
8 const Int_t    charge  = 1 ;
9 const Int_t    PDG = 2212; 
10 const Int_t    minclustersTPC = 50 ;
11 //----------------------------------------------------
12
13 Bool_t AliCFSingleTrackTaskCAF( )
14 {
15   
16   //
17   // macro to run the example task on CAF (from A.Mastroserio)
18   //
19
20   TBenchmark benchmark;
21   benchmark.Start("AliSingleTrackTask");
22   
23   //
24   // Connect to proof
25   //
26   
27   TProof::Reset("proof://user@lxb6046.cern.ch"); 
28   TProof::Open("proof://user@lxb6046.cern.ch");   
29   
30   // gProof->ClearPackage("STEERBase");
31   // gProof->ClearPackage("ESD");
32   // gProof->ClearPackage("AOD");
33   // gProof->ClearPackage("ANALYSIS");
34   
35   // Enable the STEERBase Package
36   gProof->UploadPackage("STEERBase.par");
37   gProof->EnablePackage("STEERBase");
38   // Enable the ESD Package
39   gProof->UploadPackage("ESD.par");
40   gProof->EnablePackage("ESD");
41   // Enable the AOD Package
42   gProof->UploadPackage("AOD.par");
43   gProof->EnablePackage("AOD");
44   // Enable the Analysis Package
45   gProof->UploadPackage("ANALYSIS.par");
46   gProof->EnablePackage("ANALYSIS");  
47   // Enable the CORRFW Package
48   // gProof->ClearPackage("CORRFW");
49   gProof->UploadPackage("CORRFW.par");
50   gProof->EnablePackage("CORRFW");
51   
52   gProof->ShowEnabledPackages();
53   gProof->Load("./AliCFSingleTrackTask.cxx+g");
54  
55   //
56   // Create the chain
57   //
58   gROOT->LoadMacro("CreateESDChain.C");
59   TChain* analysisChain = CreateESDChain("ESD82XX_30K.txt", 100);
60   
61   Info("AliCFSingleTrackTask",Form("CHAIN HAS %d ENTRIES",(Int_t)analysisChain->GetEntries()));
62   
63   //CONTAINER DEFINITION
64   Info("AliCFSingleTrackTask","SETUP CONTAINER");
65   //the sensitive variables (2 in this example), their indices
66   Int_t ipt = 0;
67   Int_t iy  = 1;
68   //Setting up the container grid... 
69   Int_t nstep = 4 ; //number of selection steps MC 
70   const Int_t nvar   = 2 ; //number of variables on the grid:pt,y
71   const Int_t nbin1  = 8 ; //bins in pt
72   const Int_t nbin2  = 8 ; //bins in y 
73   //arrays for the number of bins in each dimension
74   const Int_t iBin[nvar] ={nbin1,nbin2};
75   //arrays for lower bounds :
76   Double_t binLim1[nbin1+1];
77   Double_t binLim2[nbin2+1];
78   //values for bin lower bounds
79   for(Int_t i=0; i<=nbin1; i++) binLim1[i]=(Double_t)ptmin + (ptmax-ptmin)/nbin1*(Double_t)i ; 
80   for(Int_t i=0; i<=nbin2; i++) binLim2[i]=(Double_t)ymin  + (ymax-ymin)  /nbin2*(Double_t)i ;
81   //one "container" for MC
82   AliCFContainer* container = new AliCFContainer("container","container for tracks",nstep,nvar,iBin);
83   //setting the bin limits
84   container -> SetBinLimits(ipt,binLim1);
85   container -> SetBinLimits(iy,binLim2);
86
87
88   //CREATE THE  CUTS -----------------------------------------------
89   //Particle-Level cuts:  
90
91   // Gen-Level kinematic cuts
92   AliCFTrackKineCuts *mcKineCuts = new AliCFTrackKineCuts("mcKineCuts","MC-level kinematic cuts");
93   mcKineCuts->SetPtRange(ptmin,ptmax);
94   mcKineCuts->SetRapidityRange(ymin,ymax);
95   mcKineCuts->SetChargeMC(charge);
96
97   AliCFParticleGenCuts* mcGenCuts = new AliCFParticleGenCuts("mcGenCuts","MC particle generation cuts");
98   mcGenCuts->SetRequireIsPrimary();
99   mcGenCuts->SetRequirePdgCode(PDG);
100   //Acceptance Cuts
101   AliCFAcceptanceCuts *mcAccCuts = new AliCFAcceptanceCuts("mcAccCuts","MC acceptance cuts");
102   mcAccCuts->SetMinNHitITS(mintrackrefsITS);
103   mcAccCuts->SetMinNHitTPC(mintrackrefsTPC);
104
105   // Rec-Level kinematic cuts
106   AliCFTrackKineCuts *recKineCuts = new AliCFTrackKineCuts("recKineCuts","rec-level kine cuts");
107   recKineCuts->SetPtRange(ptmin,ptmax);
108   recKineCuts->SetRapidityRange(ymin,ymax);
109   recKineCuts->SetChargeRec(charge);
110   // QA histograms for rec-level kinematic cuts
111   recKineCuts->SetQAOn(kTRUE);
112
113   AliCFTrackQualityCuts *recQualityCuts = new AliCFTrackQualityCuts("recQualityCuts","rec-level quality cuts");
114   recQualityCuts->SetMinNClusterTPC(minclustersTPC);
115   recQualityCuts->SetRequireITSRefit(kTRUE);
116   // QA histograms for rec-level quality cuts
117   recQualityCuts->SetQAOn(kTRUE);
118
119   AliCFTrackIsPrimaryCuts *recIsPrimaryCuts = new AliCFTrackIsPrimaryCuts("recIsPrimaryCuts","rec-level isPrimary cuts");
120   recIsPrimaryCuts->SetMaxNSigmaToVertex(3);
121   // QA histograms for rec-level primary-check cuts
122   recIsPrimaryCuts->SetQAOn(kTRUE);
123
124   AliCFTrackCutPid* cutPID = new AliCFTrackCutPid("cutPID","ESD_PID") ;
125   Double_t prior[AliPID::kSPECIES] = {0.0244519,
126                                       0.0143988,
127                                       0.805747 ,
128                                       0.0928785,
129                                       0.0625243 };
130   cutPID->SetPriors(prior);
131   cutPID->SetProbabilityCut(0.0);
132   cutPID->SetDetectors("TPC TOF");
133   switch(TMath::Abs(PDG)) {
134   case 11   : cutPID->SetParticleType(AliPID::kElectron, kTRUE); break;
135   case 13   : cutPID->SetParticleType(AliPID::kMuon    , kTRUE); break;
136   case 211  : cutPID->SetParticleType(AliPID::kPion    , kTRUE); break;
137   case 321  : cutPID->SetParticleType(AliPID::kKaon    , kTRUE); break;
138   case 2212 : cutPID->SetParticleType(AliPID::kProton  , kTRUE); break;
139   default   : printf("UNDEFINED PID\n"); break;
140   }
141   cutPID->SetQAOn(kTRUE);
142
143   printf("CREATE MC KINE CUTS\n");
144   TObjArray* mcList = new TObjArray(0) ;
145   mcList->AddLast(mcKineCuts);
146   mcList->AddLast(mcGenCuts);
147
148   printf("CREATE ACCEPTANCE CUTS\n");
149   TObjArray* accList = new TObjArray(0) ;
150   accList->AddLast(mcAccCuts);
151
152   printf("CREATE RECONSTRUCTION CUTS\n");
153   TObjArray* recList = new TObjArray(0) ;
154   recList->AddLast(recKineCuts);
155   recList->AddLast(recQualityCuts);
156   recList->AddLast(recIsPrimaryCuts);
157
158   printf("CREATE PID CUTS\n");
159   TObjArray* fPIDCutList = new TObjArray(0) ;
160   fPIDCutList->AddLast(cutPID);
161
162
163   //CREATE THE INTERFACE TO CORRECTION FRAMEWORK USED IN THE TASK
164   printf("CREATE INTERFACE AND CUTS\n");
165   AliCFManager* man = new AliCFManager() ;
166   man->SetParticleContainer     (container);
167   man->SetParticleCutsList(AliCFManager::kPartGenCuts,mcList);
168   man->SetParticleCutsList(AliCFManager::kPartAccCuts,accList);
169   man->SetParticleCutsList(AliCFManager::kPartRecCuts,recList);
170   man->SetParticleCutsList(AliCFManager::kPartSelCuts,fPIDCutList);
171
172
173   //CREATE THE TASK
174   printf("CREATE TASK\n");
175   // create the task
176   AliCFSingleTrackTask *task = new AliCFSingleTrackTask("AliSingleTrackTask");
177   task->SetCFManager(man); //here is set the CF manager
178
179
180   //SETUP THE ANALYSIS MANAGER TO READ INPUT CHAIN AND WRITE DESIRED OUTPUTS
181   printf("CREATE ANALYSIS MANAGER\n");
182   // Make the analysis manager
183   AliAnalysisManager *mgr = new AliAnalysisManager("TestManager");
184
185   // Create and connect containers for input/output
186   //input data
187   AliAnalysisDataContainer *cinput0  = mgr->CreateContainer("cchain0",TChain::Class(),AliAnalysisManager::kInputContainer);
188   // output histo (number of events processed)
189   AliAnalysisDataContainer *coutput0 = mgr->CreateContainer("chist0", TH1I::Class(),AliAnalysisManager::kOutputContainer,"output.root");
190   // output Correction Framework Container (for acceptance & efficiency calculations)
191   AliAnalysisDataContainer *coutput1 = mgr->CreateContainer("ccontainer0", AliCFContainer::Class(),AliAnalysisManager::kOutputContainer,"container.root");
192   // output QA histograms 
193   AliAnalysisDataContainer *coutput2 = mgr->CreateContainer("clist0", TList::Class(),AliAnalysisManager::kOutputContainer,"QAhistos.root");
194
195   mgr->AddTask(task);
196   mgr->ConnectInput(task,0,cinput0);
197   mgr->ConnectOutput(task,0,coutput0);
198   mgr->ConnectOutput(task,1,coutput1);
199   mgr->ConnectOutput(task,2,coutput2);
200   cinput0->SetData(analysisChain);
201  
202   //NEW INTERFACE TO MC INFORMATION
203   AliMCEventHandler* mcHandler = new AliMCEventHandler();
204   mgr->SetMCtruthEventHandler(mcHandler);
205
206   printf("READY TO RUN\n");
207   //RUN !!!
208   if (mgr->InitAnalysis()) {
209     mgr->PrintStatus();
210     mgr->StartAnalysis("proof",analysisChain);
211   }
212
213   benchmark.Stop("AliSingleTrackTask");
214   benchmark.Show("AliSingleTrackTask");
215   
216   return kTRUE ;
217 }
218