]> git.uio.no Git - u/mrichter/AliRoot.git/blob - CORRFW/test/AliCFSingleTrackTask.C
Renaming linkdef files for consistency
[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   // SET TLIST FOR QA HISTOS
90   TList* qaList = new TList();
91
92   //CREATE THE  CUTS -----------------------------------------------
93
94   // Gen-Level kinematic cuts
95   AliCFTrackKineCuts *mcKineCuts = new AliCFTrackKineCuts("mcKineCuts","MC-level kinematic cuts");
96   mcKineCuts->SetPtRange(ptmin,ptmax);
97   mcKineCuts->SetRapidityRange(ymin,ymax);
98   mcKineCuts->SetChargeMC(charge);
99   mcKineCuts->SetQAOn(qaList);
100
101   //Particle-Level cuts:  
102   AliCFParticleGenCuts* mcGenCuts = new AliCFParticleGenCuts("mcGenCuts","MC particle generation cuts");
103   mcGenCuts->SetRequireIsPrimary();
104   mcGenCuts->SetRequirePdgCode(PDG);
105   mcGenCuts->SetQAOn(qaList);
106
107   //Acceptance Cuts
108   AliCFAcceptanceCuts *mcAccCuts = new AliCFAcceptanceCuts("mcAccCuts","MC acceptance cuts");
109   mcAccCuts->SetMinNHitITS(mintrackrefsITS);
110   mcAccCuts->SetMinNHitTPC(mintrackrefsTPC);
111   mcAccCuts->SetQAOn(qaList);
112
113   // Rec-Level kinematic cuts
114   AliCFTrackKineCuts *recKineCuts = new AliCFTrackKineCuts("recKineCuts","rec-level kine cuts");
115   recKineCuts->SetPtRange(ptmin,ptmax);
116   recKineCuts->SetRapidityRange(ymin,ymax);
117   recKineCuts->SetChargeRec(charge);
118   recKineCuts->SetQAOn(qaList);
119
120   AliCFTrackQualityCuts *recQualityCuts = new AliCFTrackQualityCuts("recQualityCuts","rec-level quality cuts");
121   recQualityCuts->SetMinNClusterTPC(minclustersTPC);
122   recQualityCuts->SetRequireITSRefit(kTRUE);
123   recQualityCuts->SetQAOn(qaList);
124
125   AliCFTrackIsPrimaryCuts *recIsPrimaryCuts = new AliCFTrackIsPrimaryCuts("recIsPrimaryCuts","rec-level isPrimary cuts");
126   recIsPrimaryCuts->SetMaxNSigmaToVertex(3);
127   recIsPrimaryCuts->SetQAOn(qaList);
128
129   AliCFTrackCutPid* cutPID = new AliCFTrackCutPid("cutPID","ESD_PID") ;
130   int n_species = AliPID::kSPECIES ;
131   Double_t* prior = new Double_t[n_species];
132   
133   prior[0] = 0.0244519 ;
134   prior[1] = 0.0143988 ;
135   prior[2] = 0.805747  ;
136   prior[3] = 0.0928785 ;
137   prior[4] = 0.0625243 ;
138   
139   cutPID->SetPriors(prior);
140   cutPID->SetProbabilityCut(0.0);
141   cutPID->SetDetectors("TPC TOF");
142   switch(TMath::Abs(PDG)) {
143   case 11   : cutPID->SetParticleType(AliPID::kElectron, kTRUE); break;
144   case 13   : cutPID->SetParticleType(AliPID::kMuon    , kTRUE); break;
145   case 211  : cutPID->SetParticleType(AliPID::kPion    , kTRUE); break;
146   case 321  : cutPID->SetParticleType(AliPID::kKaon    , kTRUE); break;
147   case 2212 : cutPID->SetParticleType(AliPID::kProton  , kTRUE); break;
148   default   : printf("UNDEFINED PID\n"); break;
149   }
150   cutPID->SetQAOn(qaList);
151
152   printf("CREATE MC KINE CUTS\n");
153   TObjArray* mcList = new TObjArray(0) ;
154   mcList->AddLast(mcKineCuts);
155   mcList->AddLast(mcGenCuts);
156
157   printf("CREATE ACCEPTANCE CUTS\n");
158   TObjArray* accList = new TObjArray(0) ;
159   accList->AddLast(mcAccCuts);
160
161   printf("CREATE RECONSTRUCTION CUTS\n");
162   TObjArray* recList = new TObjArray(0) ;
163   recList->AddLast(recKineCuts);
164   recList->AddLast(recQualityCuts);
165   recList->AddLast(recIsPrimaryCuts);
166
167   printf("CREATE PID CUTS\n");
168   TObjArray* fPIDCutList = new TObjArray(0) ;
169   fPIDCutList->AddLast(cutPID);
170
171   //CREATE THE INTERFACE TO CORRECTION FRAMEWORK USED IN THE TASK
172   printf("CREATE INTERFACE AND CUTS\n");
173   AliCFManager* man = new AliCFManager() ;
174   man->SetParticleContainer     (container);
175   man->SetParticleCutsList(AliCFManager::kPartGenCuts,mcList);
176   man->SetParticleCutsList(AliCFManager::kPartAccCuts,accList);
177   man->SetParticleCutsList(AliCFManager::kPartRecCuts,recList);
178   man->SetParticleCutsList(AliCFManager::kPartSelCuts,fPIDCutList);
179
180
181   //CREATE THE TASK
182   printf("CREATE TASK\n");
183   // create the task
184   AliCFSingleTrackTask *task = new AliCFSingleTrackTask("AliSingleTrackTask");
185   task->SetCFManager(man); //here is set the CF manager
186   task->SetQAList(qaList);
187
188   //SETUP THE ANALYSIS MANAGER TO READ INPUT CHAIN AND WRITE DESIRED OUTPUTS
189   printf("CREATE ANALYSIS MANAGER\n");
190   // Make the analysis manager
191   AliAnalysisManager *mgr = new AliAnalysisManager("TestManager");
192
193   if (useGrid) mgr->SetAnalysisType(AliAnalysisManager::kGridAnalysis);
194   else mgr->SetAnalysisType(AliAnalysisManager::kLocalAnalysis);
195
196   AliMCEventHandler*  mcHandler = new AliMCEventHandler();
197   AliESDInputHandler* esdHandler = new AliESDInputHandler();
198   mgr->SetMCtruthEventHandler(mcHandler);
199   mgr->SetInputEventHandler(esdHandler);
200
201   // Create and connect containers for input/output
202
203   //------ input data ------
204   AliAnalysisDataContainer *cinput0  = mgr->CreateContainer("cchain0",TChain::Class(),AliAnalysisManager::kInputContainer);
205
206   // ----- output data -----
207   
208   //slot 0 : default output tree (by default handled by AliAnalysisTaskSE)
209   AliAnalysisDataContainer *coutput0 = mgr->CreateContainer("ctree0", TTree::Class(),AliAnalysisManager::kOutputContainer,"output.root");
210
211   //now comes user's output objects :
212   
213   // output TH1I for event counting
214   AliAnalysisDataContainer *coutput1 = mgr->CreateContainer("chist0", TH1I::Class(),AliAnalysisManager::kOutputContainer,"output.root");
215   // output Correction Framework Container (for acceptance & efficiency calculations)
216   AliAnalysisDataContainer *coutput2 = mgr->CreateContainer("ccontainer0", AliCFContainer::Class(),AliAnalysisManager::kOutputContainer,"output.root");
217   // output QA histograms 
218   AliAnalysisDataContainer *coutput3 = mgr->CreateContainer("clist0", TList::Class(),AliAnalysisManager::kOutputContainer,"output.root");
219
220   cinput0->SetData(analysisChain);
221
222   mgr->AddTask(task);
223   mgr->ConnectInput(task,0,cinput0);
224   mgr->ConnectOutput(task,0,coutput0);
225   mgr->ConnectOutput(task,1,coutput1);
226   mgr->ConnectOutput(task,2,coutput2);
227   mgr->ConnectOutput(task,3,coutput3);
228  
229   printf("READY TO RUN\n");
230   //RUN !!!
231   if (mgr->InitAnalysis()) {
232     mgr->PrintStatus();
233     mgr->StartAnalysis("local",analysisChain);
234   }
235
236   benchmark.Stop("AliSingleTrackTask");
237   benchmark.Show("AliSingleTrackTask");
238
239   return kTRUE ;
240 }
241
242 void Load() {
243
244   //load the required aliroot libraries
245   gSystem->Load("libANALYSIS") ;
246   gSystem->Load("libANALYSISalice") ;
247   gSystem->Load("libCORRFW.so") ;
248
249   //compile online the task class
250   gSystem->SetIncludePath("-I. -I$ALICE_ROOT/include -I$ROOTSYS/include");
251   gROOT->LoadMacro("./AliCFSingleTrackTask.cxx+");
252 }