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