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