]> git.uio.no Git - u/mrichter/AliRoot.git/blob - CORRFW/test/AliCFSingleTrackTask.C
added AOD support for single tracks ; refit flags replaced by status flag in
[u/mrichter/AliRoot.git] / CORRFW / test / AliCFSingleTrackTask.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 AliCFSingleTrackTask(
14                             const Bool_t useGrid = 1,
15                             const Bool_t readAOD = 0,
16                             const Bool_t readTPCTracks = 0,
17                             const char * kTagXMLFile="wn.xml" // XML file containing tags
18                             )
19 {
20   
21   TBenchmark benchmark;
22   benchmark.Start("AliSingleTrackTask");
23
24   AliLog::SetGlobalDebugLevel(0);
25
26   Load() ; //load the required libraries
27
28   TChain * analysisChain ;
29
30   if (useGrid) { //data located on AliEn
31     TGrid::Connect("alien://") ;    //  Create an AliRunTagCuts and an AliEventTagCuts Object 
32                                     //  and impose some selection criteria
33     AliRunTagCuts      *runCuts   = new AliRunTagCuts(); 
34     AliEventTagCuts    *eventCuts = new AliEventTagCuts(); 
35     AliLHCTagCuts      *lhcCuts   = new AliLHCTagCuts(); 
36     AliDetectorTagCuts *detCuts   = new AliDetectorTagCuts(); 
37     eventCuts->SetMultiplicityRange(0,2000);
38
39     //  Create an AliTagAnalysis Object and chain the tags
40     AliTagAnalysis   *tagAna = new AliTagAnalysis(); 
41     if (readAOD) tagAna->SetType("AOD");  //for aliroot > v4-05
42     else         tagAna->SetType("ESD");  //for aliroot > v4-05
43     TAlienCollection *coll   = TAlienCollection::Open(kTagXMLFile); 
44     TGridResult      *tagResult = coll->GetGridResult("",0,0);
45     tagResult->Print();
46     tagAna->ChainGridTags(tagResult);
47
48     //  Create a new esd chain and assign the chain that is returned by querying the tags
49     analysisChain = tagAna->QueryTags(runCuts,lhcCuts,detCuts,eventCuts); 
50   }
51
52   else {// local data
53     //here put your input data path
54     printf("\n\nRunning on local file, please check the path\n\n");
55
56     if (readAOD) {
57       analysisChain = new TChain("aodTree");
58       analysisChain->Add("AliAOD.root");
59     }
60     else {
61       analysisChain = new TChain("esdTree");
62       analysisChain->Add("AliESDs.root");
63     }
64   }
65   
66
67   Info("AliCFSingleTrackTask",Form("CHAIN HAS %d ENTRIES",(Int_t)analysisChain->GetEntries()));
68
69   //CONTAINER DEFINITION
70   Info("AliCFSingleTrackTask","SETUP CONTAINER");
71   //the sensitive variables (2 in this example), their indices
72   UInt_t ipt = 0;
73   UInt_t iy  = 1;
74   //Setting up the container grid... 
75   UInt_t nstep = 4 ; //number of selection steps MC 
76   const Int_t nvar   = 2 ; //number of variables on the grid:pt,y
77   const Int_t nbin1  = 8 ; //bins in pt
78   const Int_t nbin2  = 8 ; //bins in y 
79
80   //arrays for the number of bins in each dimension
81   Int_t iBin[nvar];
82   iBin[0]=nbin1;
83   iBin[1]=nbin2;
84
85   //arrays for lower bounds :
86   Double_t *binLim1=new Double_t[nbin1+1];
87   Double_t *binLim2=new Double_t[nbin2+1];
88
89   //values for bin lower bounds
90   for(Int_t i=0; i<=nbin1; i++) binLim1[i]=(Double_t)ptmin + (ptmax-ptmin)/nbin1*(Double_t)i ; 
91   for(Int_t i=0; i<=nbin2; i++) binLim2[i]=(Double_t)ymin  + (ymax-ymin)  /nbin2*(Double_t)i ;
92   //one "container" for MC
93   AliCFContainer* container = new AliCFContainer("container","container for tracks",nstep,nvar,iBin);
94   //setting the bin limits
95   container -> SetBinLimits(ipt,binLim1);
96   container -> SetBinLimits(iy,binLim2);
97
98
99   // SET TLIST FOR QA HISTOS
100   TList* qaList = new TList();
101
102   //CREATE THE  CUTS -----------------------------------------------
103
104   // Gen-Level kinematic cuts
105   AliCFTrackKineCuts *mcKineCuts = new AliCFTrackKineCuts("mcKineCuts","MC-level kinematic cuts");
106   mcKineCuts->SetPtRange(ptmin,ptmax);
107   mcKineCuts->SetRapidityRange(ymin,ymax);
108   mcKineCuts->SetChargeMC(charge);
109   mcKineCuts->SetQAOn(qaList);
110
111   //Particle-Level cuts:  
112   AliCFParticleGenCuts* mcGenCuts = new AliCFParticleGenCuts("mcGenCuts","MC particle generation cuts");
113   mcGenCuts->SetRequireIsPrimary();
114   mcGenCuts->SetRequirePdgCode(PDG);
115   mcGenCuts->SetQAOn(qaList);
116
117   //Acceptance Cuts
118   AliCFAcceptanceCuts *mcAccCuts = new AliCFAcceptanceCuts("mcAccCuts","MC acceptance cuts");
119   mcAccCuts->SetMinNHitITS(mintrackrefsITS);
120   mcAccCuts->SetMinNHitTPC(mintrackrefsTPC);
121   mcAccCuts->SetQAOn(qaList);
122
123   // Rec-Level kinematic cuts
124   AliCFTrackKineCuts *recKineCuts = new AliCFTrackKineCuts("recKineCuts","rec-level kine cuts");
125   recKineCuts->SetPtRange(ptmin,ptmax);
126   recKineCuts->SetRapidityRange(ymin,ymax);
127   recKineCuts->SetChargeRec(charge);
128   recKineCuts->SetQAOn(qaList);
129
130   AliCFTrackQualityCuts *recQualityCuts = new AliCFTrackQualityCuts("recQualityCuts","rec-level quality cuts");
131   if (!readAOD)       recQualityCuts->SetMinNClusterTPC(minclustersTPC);
132   if (!readTPCTracks) recQualityCuts->SetStatus(AliESDtrack::kITSrefit);
133   recQualityCuts->SetQAOn(qaList);
134
135   AliCFTrackIsPrimaryCuts *recIsPrimaryCuts = new AliCFTrackIsPrimaryCuts("recIsPrimaryCuts","rec-level isPrimary cuts");
136   if (readAOD) recIsPrimaryCuts->SetAODType(AliAODTrack::kPrimary);
137   else         recIsPrimaryCuts->SetMaxNSigmaToVertex(3);
138   recIsPrimaryCuts->SetQAOn(qaList);
139
140   AliCFTrackCutPid* cutPID = new AliCFTrackCutPid("cutPID","ESD_PID") ;
141   int n_species = AliPID::kSPECIES ;
142   Double_t* prior = new Double_t[n_species];
143   
144   prior[0] = 0.0244519 ;
145   prior[1] = 0.0143988 ;
146   prior[2] = 0.805747  ;
147   prior[3] = 0.0928785 ;
148   prior[4] = 0.0625243 ;
149   
150   cutPID->SetPriors(prior);
151   cutPID->SetProbabilityCut(0.0);
152   if (readTPCTracks) cutPID->SetDetectors("TPC");
153   else               cutPID->SetDetectors("ALL");
154   if (readAOD) cutPID->SetAODmode(kTRUE );
155   else         cutPID->SetAODmode(kFALSE);
156   switch(TMath::Abs(PDG)) {
157   case 11   : cutPID->SetParticleType(AliPID::kElectron, kTRUE); break;
158   case 13   : cutPID->SetParticleType(AliPID::kMuon    , kTRUE); break;
159   case 211  : cutPID->SetParticleType(AliPID::kPion    , kTRUE); break;
160   case 321  : cutPID->SetParticleType(AliPID::kKaon    , kTRUE); break;
161   case 2212 : cutPID->SetParticleType(AliPID::kProton  , kTRUE); break;
162   default   : printf("UNDEFINED PID\n"); break;
163   }
164   cutPID->SetQAOn(qaList);
165
166   printf("CREATE MC KINE CUTS\n");
167   TObjArray* mcList = new TObjArray(0) ;
168   mcList->AddLast(mcKineCuts);
169   mcList->AddLast(mcGenCuts);
170
171   printf("CREATE ACCEPTANCE CUTS\n");
172   TObjArray* accList = new TObjArray(0) ;
173   accList->AddLast(mcAccCuts);
174
175   printf("CREATE RECONSTRUCTION CUTS\n");
176   TObjArray* recList = new TObjArray(0) ;
177   recList->AddLast(recKineCuts);
178   recList->AddLast(recQualityCuts);
179   recList->AddLast(recIsPrimaryCuts);
180
181   printf("CREATE PID CUTS\n");
182   TObjArray* fPIDCutList = new TObjArray(0) ;
183   fPIDCutList->AddLast(cutPID);
184
185   //CREATE THE INTERFACE TO CORRECTION FRAMEWORK USED IN THE TASK
186   printf("CREATE INTERFACE AND CUTS\n");
187   AliCFManager* man = new AliCFManager() ;
188   man->SetParticleContainer     (container);
189   man->SetParticleCutsList(AliCFManager::kPartGenCuts,mcList);
190   man->SetParticleCutsList(AliCFManager::kPartAccCuts,accList);
191   man->SetParticleCutsList(AliCFManager::kPartRecCuts,recList);
192   man->SetParticleCutsList(AliCFManager::kPartSelCuts,fPIDCutList);
193
194
195   //CREATE THE TASK
196   printf("CREATE TASK\n");
197   // create the task
198   AliCFSingleTrackTask *task = new AliCFSingleTrackTask("AliSingleTrackTask");
199   task->SetCFManager(man); //here is set the CF manager
200   task->SetQAList(qaList);
201   if (readAOD)       task->SetReadAODData() ;
202   if (readTPCTracks) task->SetReadTPCTracks();
203
204   //SETUP THE ANALYSIS MANAGER TO READ INPUT CHAIN AND WRITE DESIRED OUTPUTS
205   printf("CREATE ANALYSIS MANAGER\n");
206   // Make the analysis manager
207   AliAnalysisManager *mgr = new AliAnalysisManager("TestManager");
208
209   if (useGrid) mgr->SetAnalysisType(AliAnalysisManager::kGridAnalysis);
210   else mgr->SetAnalysisType(AliAnalysisManager::kLocalAnalysis);
211
212
213   AliMCEventHandler*  mcHandler = new AliMCEventHandler();
214   mgr->SetMCtruthEventHandler(mcHandler);
215  
216   AliInputEventHandler* dataHandler ;
217   
218   if   (readAOD) dataHandler = new AliAODInputHandler();
219   else           dataHandler = new AliESDInputHandler();
220   mgr->SetInputEventHandler(dataHandler);
221
222   // Create and connect containers for input/output
223
224   //------ input data ------
225   AliAnalysisDataContainer *cinput0  = mgr->CreateContainer("cchain0",TChain::Class(),AliAnalysisManager::kInputContainer);
226
227   // ----- output data -----
228   
229   //slot 0 : default output tree (by default handled by AliAnalysisTaskSE)
230   AliAnalysisDataContainer *coutput0 = mgr->CreateContainer("ctree0", TTree::Class(),AliAnalysisManager::kOutputContainer,"output.root");
231
232   //now comes user's output objects :
233   
234   // output TH1I for event counting
235   AliAnalysisDataContainer *coutput1 = mgr->CreateContainer("chist0", TH1I::Class(),AliAnalysisManager::kOutputContainer,"output.root");
236   // output Correction Framework Container (for acceptance & efficiency calculations)
237   AliAnalysisDataContainer *coutput2 = mgr->CreateContainer("ccontainer0", AliCFContainer::Class(),AliAnalysisManager::kOutputContainer,"output.root");
238   // output QA histograms 
239   AliAnalysisDataContainer *coutput3 = mgr->CreateContainer("clist0", TList::Class(),AliAnalysisManager::kOutputContainer,"output.root");
240
241   cinput0->SetData(analysisChain);
242
243   mgr->AddTask(task);
244   mgr->ConnectInput(task,0,cinput0);
245   mgr->ConnectOutput(task,0,coutput0);
246   mgr->ConnectOutput(task,1,coutput1);
247   mgr->ConnectOutput(task,2,coutput2);
248   mgr->ConnectOutput(task,3,coutput3);
249  
250   printf("READY TO RUN\n");
251   //RUN !!!
252   if (mgr->InitAnalysis()) {
253     mgr->PrintStatus();
254     mgr->StartAnalysis("local",analysisChain);
255   }
256
257   benchmark.Stop("AliSingleTrackTask");
258   benchmark.Show("AliSingleTrackTask");
259
260   return kTRUE ;
261 }
262
263 void Load() {
264
265   //load the required aliroot libraries
266   gSystem->Load("libANALYSIS") ;
267   gSystem->Load("libANALYSISalice") ;
268   gSystem->Load("libCORRFW.so") ;
269
270   //compile online the task class
271   gSystem->SetIncludePath("-I. -I$ALICE_ROOT/include -I$ALICE_ROOT/TPC -I$ROOTSYS/include");
272   gROOT->LoadMacro("./AliCFSingleTrackTask.cxx+");
273 }