]> git.uio.no Git - u/mrichter/AliRoot.git/blob - CORRFW/test/AliCFV0Task.C
Updated branch aliroot-master, development and master to
[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("your_data_path/001/AliAOD.root");
56       analysisChain->Add("your_data_path/002/AliAOD.root");
57     }
58     else {
59       analysisChain = new TChain("esdTree");
60       analysisChain->Add("your_data_path/001/AliESDs.root");
61       analysisChain->Add("your_data_path/002/AliESDs.root");
62     }
63   }
64   
65   
66   Info("AliCFV0Task",Form("CHAIN HAS %d ENTRIES",(Int_t)analysisChain->GetEntries()));
67
68   //CONTAINER DEFINITION
69   Info("AliCFV0Task","SETUP CONTAINER");
70   //the sensitive variables, their indices
71   Int_t ipt = 0;
72   Int_t iy  = 1;
73   //Setting up the container grid... 
74   Int_t nstep = 4 ; //number of selection steps MC 
75   const Int_t nvar   = 2 ; //number of variables on the grid:pt,y,phi,vtx
76   const Int_t nbin1  = 8 ; //bins in pt
77   const Int_t nbin2  = 8 ; //bins in y 
78   //arrays for the number of bins in each dimension
79   const Int_t iBin[nvar] ={nbin1,nbin2};
80   //arrays for lower bounds :
81   Double_t binLim1[nbin1+1];
82   Double_t binLim2[nbin2+1];
83   //values for bin lower bounds
84   for(Int_t i=0; i<=nbin1; i++) binLim1[i]=ptmin + (ptmax-ptmin)/nbin1*(Double_t)i ; 
85   for(Int_t i=0; i<=nbin2; i++) binLim2[i]=ymin  + (ymax-ymin)  /nbin2*(Double_t)i ;
86   //one "container" for MC
87   AliCFContainer* container = new AliCFContainer("container","container for V0s",nstep,nvar,iBin);
88   //setting the bin limits
89   container -> SetBinLimits(ipt,binLim1);
90   container -> SetBinLimits(iy,binLim2);
91
92
93   //CREATE THE  CUTS -----------------------------------------------
94   
95   // Gen-Level kinematic cuts
96   AliCFTrackKineCuts *mcKineCuts = new AliCFTrackKineCuts("mcKineCuts","MC-level kinematic cuts");
97   mcKineCuts->SetPtRange(ptmin,ptmax);
98   mcKineCuts->SetRapidityRange(ymin,ymax);
99   mcKineCuts->SetChargeMC(0);
100
101   AliCFParticleGenCuts* mcGenCuts = new AliCFParticleGenCuts("mcGenCuts","MC particle generation cuts");
102   //mcGenCuts->SetRequireIsPrimary(); //problem with some particles...
103   mcGenCuts->SetRequirePdgCode(PDG);
104
105   //Acceptance Cuts
106   AliCFPairAcceptanceCuts *mcAccCuts = new AliCFPairAcceptanceCuts("mcAccCuts","MC acceptance cuts for V0");
107   mcAccCuts->GetNegCut()->SetMinNHitITS(mintrackrefsITS);
108   mcAccCuts->GetPosCut()->SetMinNHitITS(mintrackrefsITS);
109   mcAccCuts->GetNegCut()->SetMinNHitTPC(mintrackrefsTPC);
110   mcAccCuts->GetPosCut()->SetMinNHitTPC(mintrackrefsTPC);
111
112   // Rec-Level kinematic cuts
113   AliCFTrackKineCuts *recKineCuts = new AliCFTrackKineCuts("recKineCuts","V0 rec-level kine cuts");
114   recKineCuts->SetPtRange(ptmin,ptmax);
115   recKineCuts->SetRapidityRange(ymin,ymax);
116   recKineCuts->SetChargeRec(0);
117
118   AliCFPairQualityCuts *recQualityCuts = new AliCFPairQualityCuts("recQualityCuts","V0 rec-level quality cuts");
119   if (!readAOD) {
120     recQualityCuts->GetNegCut()->SetMinNClusterTPC(minclustersTPC);
121     recQualityCuts->GetPosCut()->SetMinNClusterTPC(minclustersTPC);
122   }
123   recQualityCuts->GetNegCut()->SetStatus(AliESDtrack::kTPCrefit & AliESDtrack::kITSrefit);
124   recQualityCuts->GetPosCut()->SetStatus(AliESDtrack::kTPCrefit & AliESDtrack::kITSrefit);
125
126   AliCFV0TopoCuts *recTopoCuts = new AliCFV0TopoCuts("recTopoCuts","V0 Topological Cuts");
127   recTopoCuts->SetMaxDcaDaughters(0.1);
128   recTopoCuts->SetMinCosPointAngle(0.995);
129   recTopoCuts->SetMinDcaNeg(0.1);
130   recTopoCuts->SetMinDcaPos(0.1);
131
132
133   AliCFPairPidCut* cutPID = new AliCFPairPidCut("cutPID","ESD_PID") ;
134   Double_t prior_pp[AliPID::kSPECIES] = {0.0244519,
135                                          0.0143988,
136                                          0.805747 ,
137                                          0.0928785,
138                                          0.0625243 };
139
140   Double_t prior_pbpb[AliPID::kSPECIES] = {0.0609,
141                                            0.1064,
142                                            0.7152 ,
143                                            0.0442,
144                                            0.0733 };
145   cutPID->GetNegCut()->SetPriors(prior_pp);
146   cutPID->GetPosCut()->SetPriors(prior_pp);
147   cutPID->GetNegCut()->SetDetectors("ITS TPC TRD");
148   cutPID->GetPosCut()->SetDetectors("ITS TPC TRD");
149
150   if (readAOD) {
151     cutPID->GetNegCut()->SetAODmode(kTRUE);
152     cutPID->GetPosCut()->SetAODmode(kTRUE);
153   }
154   else {
155     cutPID->GetNegCut()->SetAODmode(kFALSE);
156     cutPID->GetPosCut()->SetAODmode(kFALSE);
157   }
158
159   cutPID->GetNegCut()->SetProbabilityCut(0);
160   cutPID->GetPosCut()->SetProbabilityCut(0);
161
162   switch(TMath::Abs(PDG)) {
163   case  310    : 
164     cutPID->GetNegCut()->SetParticleType(AliPID::kPion  , kTRUE);
165     cutPID->GetPosCut()->SetParticleType(AliPID::kPion  , kTRUE);
166     break;
167   case  3122   : 
168     cutPID->GetNegCut()->SetParticleType(AliPID::kPion  , kTRUE);
169     cutPID->GetPosCut()->SetParticleType(AliPID::kProton, kTRUE);
170     break;
171   case -3122   : 
172     cutPID->GetNegCut()->SetParticleType(AliPID::kProton, kTRUE);
173     cutPID->GetPosCut()->SetParticleType(AliPID::kPion  , kTRUE);
174     break;
175   default      : printf("UNDEFINED PID\n"); break;
176   }
177
178   Info("AliCFV0Task","CREATE MC KINE CUTS");
179   TObjArray* mcList = new TObjArray(0) ;
180   mcList->AddLast(mcKineCuts);
181   mcList->AddLast(mcGenCuts);
182
183   Info("AliCFV0Task","CREATE ACCEPTANCE CUTS");
184   TObjArray* accList = new TObjArray(0) ;
185   accList->AddLast(mcAccCuts);
186
187   Info("AliCFV0Task","CREATE RECONSTRUCTION CUTS");
188   TObjArray* recList = new TObjArray(0) ;
189   recList->AddLast(recKineCuts);
190   recList->AddLast(recQualityCuts);
191   recList->AddLast(recTopoCuts);
192
193   Info("AliCFV0Task","CREATE PID CUTS");
194   TObjArray* fPIDCutList = new TObjArray(0) ;
195   fPIDCutList->AddLast(cutPID);
196
197
198   //CREATE THE INTERFACE TO CORRECTION FRAMEWORK USED IN THE TASK
199   Info("AliCFV0Task","CREATE INTERFACE AND CUTS");
200   AliCFManager* man = new AliCFManager() ;
201   man->SetParticleContainer(container);
202   man->SetParticleCutsList (AliCFManager::kPartGenCuts,mcList);
203   man->SetParticleCutsList (AliCFManager::kPartAccCuts,accList);
204   man->SetParticleCutsList (AliCFManager::kPartRecCuts,recList);
205   man->SetParticleCutsList (AliCFManager::kPartSelCuts,fPIDCutList);
206
207   //CREATE THE TASK
208   Info("AliCFV0Task","CREATE TASK");
209   // create the task
210   AliCFV0Task *task = new AliCFV0Task("AliCFV0Task");
211   task->SetCFManager(man); //here is set the CF manager
212   if (!readAOD)      task->SetRebuildV0s(kTRUE);
213   task->SetV0PDG(PDG);
214
215   //SETUP THE ANALYSIS MANAGER TO READ INPUT CHAIN AND WRITE DESIRED OUTPUTS
216   Info("AliCFV0Task","CREATE ANALYSIS MANAGER");
217   // Make the analysis manager
218   AliAnalysisManager *mgr = new AliAnalysisManager("TestManager");
219
220   if (useGrid) mgr->SetAnalysisType(AliAnalysisManager::kGridAnalysis);
221   else mgr->SetAnalysisType(AliAnalysisManager::kLocalAnalysis);
222
223   AliMCEventHandler* mcHandler = new AliMCEventHandler();
224   mgr->SetMCtruthEventHandler(mcHandler);
225
226   AliInputEventHandler* dataHandler ;
227
228   if   (readAOD) dataHandler = new AliAODInputHandler();
229   else           dataHandler = new AliESDInputHandler();
230   mgr->SetInputEventHandler(dataHandler);
231
232   // Create and connect containers for input/output
233
234   //input data
235   AliAnalysisDataContainer *cinput0  = 
236     mgr->CreateContainer("cchain0",TChain::Class(),AliAnalysisManager::kInputContainer);
237   //slot 0 : default output tree (by default handled by AliAnalysisTaskSE)
238   AliAnalysisDataContainer *coutput0 = mgr->CreateContainer("ctree0", TTree::Class(),AliAnalysisManager::kOutputContainer,"outputV0.root");
239   // output histo (number of events processed)
240   AliAnalysisDataContainer *coutput1 = 
241     mgr->CreateContainer("chist0", TH1I::Class(),AliAnalysisManager::kOutputContainer,"outputV0.root");
242   // output Correction Framework Container (for acceptance & efficiency calculations)
243   AliAnalysisDataContainer *coutput2 = 
244     mgr->CreateContainer("ccontainer0", AliCFContainer::Class(),AliAnalysisManager::kOutputContainer,"outputV0.root");
245
246   cinput0->SetData(analysisChain);
247
248   mgr->AddTask(task);
249   mgr->ConnectInput (task,0,mgr->GetCommonInputContainer());
250   mgr->ConnectOutput(task,0,coutput0);
251   mgr->ConnectOutput(task,1,coutput1);
252   mgr->ConnectOutput(task,2,coutput2);
253
254  
255   Info("AliCFV0Task","READY TO RUN");
256   //RUN !!!
257   if (mgr->InitAnalysis()) {
258     mgr->PrintStatus();
259     mgr->StartAnalysis("local",analysisChain);
260   }
261
262   benchmark.Stop("AliCFV0Task");
263   benchmark.Show("AliCFV0Task");
264   
265   return kTRUE ;
266 }
267
268 void Load(Bool_t useGrid) {
269
270   //load the required aliroot libraries
271   gSystem->Load("libANALYSIS") ;
272   gSystem->Load("libANALYSISalice") ;
273   gSystem->Load("libCORRFW.so") ;
274
275   //compile online the task class
276   gSystem->SetIncludePath("-I. -I$ALICE_ROOT/include -I$ROOTSYS/include");
277   gROOT->LoadMacro("./AliCFV0Task.cxx+");
278 }