]> git.uio.no Git - u/mrichter/AliRoot.git/blob - CORRFW/test/AliCFV0Task.C
f8556f2c31e89a010fcc760238d3074da54ceafe
[u/mrichter/AliRoot.git] / CORRFW / test / AliCFV0Task.C
1 //DEFINITION OF A FEW CONSTANTS
2 const Double_t ymin  = -2.0 ;
3 const Double_t ymax  =  2.0 ;
4 const Double_t ptmin =  0.0 ;
5 const Double_t ptmax =  8.0 ;
6 const Double_t phimin =   0 ;
7 const Double_t phimax = 360 ;
8 const Int_t    mintrackrefsTPC = 2 ;
9 const Int_t    mintrackrefsITS = 3 ;
10 const Int_t    PDG = 310; 
11 const Int_t    minclustersTPC = 10 ;
12 const Int_t    chargeV0 = 0 ;
13 //----------------------------------------------------
14
15 Bool_t AliCFV0Task(
16                    const Bool_t useGrid = 1,
17                    const Bool_t readAOD = 0,
18                    const char * kTagXMLFile="wn.xml", // XML file containing tags
19                    )
20 {
21   
22   TBenchmark benchmark;
23   benchmark.Start("AliSingleTrackTask");
24
25   AliLog::SetGlobalDebugLevel(0);
26
27   Load(useGrid) ; //load the required libraries
28
29   AliLog::SetGlobalDebugLevel(0);
30
31   TChain * analysisChain ;
32
33   if (useGrid) { //data located on AliEn
34     TGrid::Connect("alien://") ;
35     //  Create an AliRunTagCuts and an AliEventTagCuts Object and impose some selection criteria
36     AliRunTagCuts      *runCuts   = new AliRunTagCuts(); 
37     AliEventTagCuts    *eventCuts = new AliEventTagCuts(); 
38     AliLHCTagCuts      *lhcCuts   = new AliLHCTagCuts(); 
39     AliDetectorTagCuts *detCuts   = new AliDetectorTagCuts(); 
40     eventCuts->SetMultiplicityRange(0,20000);
41     //  Create an AliTagAnalysis Object and chain the tags
42     AliTagAnalysis   *tagAna = new AliTagAnalysis(); 
43     tagAna->SetType("ESD");  //for aliroot > v4-05
44     TAlienCollection *coll   = TAlienCollection::Open(kTagXMLFile); 
45     TGridResult      *tagResult = coll->GetGridResult("",0,0);
46     tagResult->Print();
47     tagAna->ChainGridTags(tagResult);
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   else {// local data
52     //here put your input data path
53     if (readAOD) {
54       analysisChain = new TChain("aodTree");
55       analysisChain->Add("AliAOD.root");
56     }
57     else {
58       analysisChain = new TChain("esdTree");
59       analysisChain->Add("AliESDs.root");
60     }
61   }
62   
63   
64   Info("AliCFV0Task",Form("CHAIN HAS %d ENTRIES",(Int_t)analysisChain->GetEntries()));
65
66   //CONTAINER DEFINITION
67   Info("AliCFV0Task","SETUP CONTAINER");
68   //the sensitive variables, their indices
69   Int_t ipt = 0;
70   Int_t iy  = 1;
71   //Setting up the container grid... 
72   Int_t nstep = 4 ; //number of selection steps MC 
73   const Int_t nvar   = 2 ; //number of variables on the grid:pt,y,phi,vtx
74   const Int_t nbin1  = 8 ; //bins in pt
75   const Int_t nbin2  = 8 ; //bins in y 
76   //arrays for the number of bins in each dimension
77   const Int_t iBin[nvar] ={nbin1,nbin2};
78   //arrays for lower bounds :
79   Double_t binLim1[nbin1+1];
80   Double_t binLim2[nbin2+1];
81   //values for bin lower bounds
82   for(Int_t i=0; i<=nbin1; i++) binLim1[i]=ptmin + (ptmax-ptmin)/nbin1*(Double_t)i ; 
83   for(Int_t i=0; i<=nbin2; i++) binLim2[i]=ymin  + (ymax-ymin)  /nbin2*(Double_t)i ;
84   //one "container" for MC
85   AliCFContainer* container = new AliCFContainer("container","container for V0s",nstep,nvar,iBin);
86   //setting the bin limits
87   container -> SetBinLimits(ipt,binLim1);
88   container -> SetBinLimits(iy,binLim2);
89
90
91   //CREATE THE  CUTS -----------------------------------------------
92   
93   // Gen-Level kinematic cuts
94   AliCFTrackKineCuts *mcKineCuts = new AliCFTrackKineCuts("mcKineCuts","MC-level kinematic cuts");
95   mcKineCuts->SetPtRange(ptmin,ptmax);
96   mcKineCuts->SetRapidityRange(ymin,ymax);
97   mcKineCuts->SetChargeMC(0);
98
99   AliCFParticleGenCuts* mcGenCuts = new AliCFParticleGenCuts("mcGenCuts","MC particle generation cuts");
100   //mcGenCuts->SetRequireIsPrimary(); //problem with some particles...
101   mcGenCuts->SetRequirePdgCode(PDG);
102
103   //Acceptance Cuts
104   AliCFPairAcceptanceCuts *mcAccCuts = new AliCFPairAcceptanceCuts("mcAccCuts","MC acceptance cuts for V0");
105   mcAccCuts->SetMinNHitITS(mintrackrefsITS,mintrackrefsITS);
106   mcAccCuts->SetMinNHitTPC(mintrackrefsTPC,mintrackrefsTPC);
107
108   // Rec-Level kinematic cuts
109   AliCFTrackKineCuts *recKineCuts = new AliCFTrackKineCuts("recKineCuts","V0 rec-level kine cuts");
110   recKineCuts->SetPtRange(ptmin,ptmax);
111   recKineCuts->SetRapidityRange(ymin,ymax);
112   recKineCuts->SetChargeRec(0);
113
114   AliCFPairQualityCuts *recQualityCuts = new AliCFPairQualityCuts("recQualityCuts","V0 rec-level quality cuts");
115   if (!readAOD) recQualityCuts->SetMinNClusterTPC(minclustersTPC,minclustersTPC);
116   recQualityCuts->SetStatus(AliESDtrack::kTPCrefit & AliESDtrack::kITSrefit, 
117                             AliESDtrack::kTPCrefit & AliESDtrack::kITSrefit) ;
118
119   AliCFV0TopoCuts *recTopoCuts = new AliCFV0TopoCuts("recTopoCuts","V0 Topological Cuts");
120   recTopoCuts->SetMaxDcaDaughters(0.1);
121   recTopoCuts->SetMinCosPointAngle(0.995);
122   recTopoCuts->SetMinDcaNeg(0.1);
123   recTopoCuts->SetMinDcaPos(0.1);
124
125
126   AliCFPairPidCut* cutPID = new AliCFPairPidCut("cutPID","ESD_PID") ;
127   Double_t prior_pp[AliPID::kSPECIES] = {0.0244519,
128                                          0.0143988,
129                                          0.805747 ,
130                                          0.0928785,
131                                          0.0625243 };
132
133   Double_t prior_pbpb[AliPID::kSPECIES] = {0.0609,
134                                            0.1064,
135                                            0.7152 ,
136                                            0.0442,
137                                            0.0733 };
138   cutPID->SetPriors(prior_pp);
139   cutPID->SetDetectors("ITS TPC TRD","ITS TPC TRD");
140   if (readAOD) cutPID->SetAODmode(kTRUE);
141   else         cutPID->SetAODmode(kFALSE);
142   cutPID->SetProbabilityCut(0,0);
143   switch(TMath::Abs(PDG)) {
144   case  310    : cutPID->SetParticleType(AliPID::kPion  , kTRUE, AliPID::kPion  , kTRUE); break;
145   case  3122   : cutPID->SetParticleType(AliPID::kPion  , kTRUE, AliPID::kProton, kTRUE); break;
146   case -3122   : cutPID->SetParticleType(AliPID::kProton, kTRUE, AliPID::kPion  , kTRUE); break;
147   default      : printf("UNDEFINED PID\n"); break;
148   }
149
150   Info("AliCFV0Task","CREATE MC KINE CUTS");
151   TObjArray* mcList = new TObjArray(0) ;
152   mcList->AddLast(mcKineCuts);
153   mcList->AddLast(mcGenCuts);
154
155   Info("AliCFV0Task","CREATE ACCEPTANCE CUTS");
156   TObjArray* accList = new TObjArray(0) ;
157   accList->AddLast(mcAccCuts);
158
159   Info("AliCFV0Task","CREATE RECONSTRUCTION CUTS");
160   TObjArray* recList = new TObjArray(0) ;
161   recList->AddLast(recKineCuts);
162   recList->AddLast(recQualityCuts);
163   recList->AddLast(recTopoCuts);
164
165   Info("AliCFV0Task","CREATE PID CUTS");
166   TObjArray* fPIDCutList = new TObjArray(0) ;
167   fPIDCutList->AddLast(cutPID);
168
169
170   //CREATE THE INTERFACE TO CORRECTION FRAMEWORK USED IN THE TASK
171   Info("AliCFV0Task","CREATE INTERFACE AND CUTS");
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   //CREATE THE TASK
180   Info("AliCFV0Task","CREATE TASK");
181   // create the task
182   AliCFV0Task *task = new AliCFV0Task("AliCFV0Task");
183   task->SetCFManager(man); //here is set the CF manager
184   if (!readAOD)      task->SetRebuildV0s(kTRUE);
185   task->SetV0PDG(PDG);
186
187   //SETUP THE ANALYSIS MANAGER TO READ INPUT CHAIN AND WRITE DESIRED OUTPUTS
188   Info("AliCFV0Task","CREATE ANALYSIS MANAGER");
189   // Make the analysis manager
190   AliAnalysisManager *mgr = new AliAnalysisManager("TestManager");
191
192   if (useGrid) mgr->SetAnalysisType(AliAnalysisManager::kGridAnalysis);
193   else mgr->SetAnalysisType(AliAnalysisManager::kLocalAnalysis);
194
195   AliMCEventHandler* mcHandler = new AliMCEventHandler();
196   mgr->SetMCtruthEventHandler(mcHandler);
197
198   AliInputEventHandler* dataHandler ;
199
200   if   (readAOD) dataHandler = new AliAODInputHandler();
201   else           dataHandler = new AliESDInputHandler();
202   mgr->SetInputEventHandler(dataHandler);
203
204   // Create and connect containers for input/output
205
206   //input data
207   AliAnalysisDataContainer *cinput0  = 
208     mgr->CreateContainer("cchain0",TChain::Class(),AliAnalysisManager::kInputContainer);
209   //slot 0 : default output tree (by default handled by AliAnalysisTaskSE)
210   AliAnalysisDataContainer *coutput0 = mgr->CreateContainer("ctree0", TTree::Class(),AliAnalysisManager::kOutputContainer,"outputV0.root");
211   // output histo (number of events processed)
212   AliAnalysisDataContainer *coutput1 = 
213     mgr->CreateContainer("chist0", TH1I::Class(),AliAnalysisManager::kOutputContainer,"outputV0.root");
214   // output Correction Framework Container (for acceptance & efficiency calculations)
215   AliAnalysisDataContainer *coutput2 = 
216     mgr->CreateContainer("ccontainer0", AliCFContainer::Class(),AliAnalysisManager::kOutputContainer,"outputV0.root");
217
218   cinput0->SetData(analysisChain);
219
220   mgr->AddTask(task);
221   mgr->ConnectInput (task,0,mgr->GetCommonInputContainer());
222   mgr->ConnectOutput(task,0,coutput0);
223   mgr->ConnectOutput(task,1,coutput1);
224   mgr->ConnectOutput(task,2,coutput2);
225
226  
227   Info("AliCFV0Task","READY TO RUN");
228   //RUN !!!
229   if (mgr->InitAnalysis()) {
230     mgr->PrintStatus();
231     mgr->StartAnalysis("local",analysisChain);
232   }
233
234   benchmark.Stop("AliCFV0Task");
235   benchmark.Show("AliCFV0Task");
236   
237   return kTRUE ;
238 }
239
240 void Load(Bool_t useGrid) {
241
242   //load the required aliroot libraries
243   gSystem->Load("libANALYSIS") ;
244   gSystem->Load("libANALYSISalice") ;
245   gSystem->Load("libCORRFW.so") ;
246
247   //compile online the task class
248   gSystem->SetIncludePath("-I. -I$ALICE_ROOT/include -I$ROOTSYS/include");
249   gROOT->LoadMacro("./AliCFV0Task.cxx+");
250 }