]> git.uio.no Git - u/mrichter/AliRoot.git/blob - CORRFW/test/AliCFV0Task.C
added a few comments
[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->GetNegCut()->SetMinNHitITS(mintrackrefsITS);
106   mcAccCuts->GetPosCut()->SetMinNHitITS(mintrackrefsITS);
107   mcAccCuts->GetNegCut()->SetMinNHitTPC(mintrackrefsTPC);
108   mcAccCuts->GetPosCut()->SetMinNHitTPC(mintrackrefsTPC);
109
110   // Rec-Level kinematic cuts
111   AliCFTrackKineCuts *recKineCuts = new AliCFTrackKineCuts("recKineCuts","V0 rec-level kine cuts");
112   recKineCuts->SetPtRange(ptmin,ptmax);
113   recKineCuts->SetRapidityRange(ymin,ymax);
114   recKineCuts->SetChargeRec(0);
115
116   AliCFPairQualityCuts *recQualityCuts = new AliCFPairQualityCuts("recQualityCuts","V0 rec-level quality cuts");
117   if (!readAOD) {
118     recQualityCuts->GetNegCut()->SetMinNClusterTPC(minclustersTPC);
119     recQualityCuts->GetPosCut()->SetMinNClusterTPC(minclustersTPC);
120   }
121   recQualityCuts->GetNegCut()->SetStatus(AliESDtrack::kTPCrefit & AliESDtrack::kITSrefit);
122   recQualityCuts->GetPosCut()->SetStatus(AliESDtrack::kTPCrefit & AliESDtrack::kITSrefit);
123
124   AliCFV0TopoCuts *recTopoCuts = new AliCFV0TopoCuts("recTopoCuts","V0 Topological Cuts");
125   recTopoCuts->SetMaxDcaDaughters(0.1);
126   recTopoCuts->SetMinCosPointAngle(0.995);
127   recTopoCuts->SetMinDcaNeg(0.1);
128   recTopoCuts->SetMinDcaPos(0.1);
129
130
131   AliCFPairPidCut* cutPID = new AliCFPairPidCut("cutPID","ESD_PID") ;
132   Double_t prior_pp[AliPID::kSPECIES] = {0.0244519,
133                                          0.0143988,
134                                          0.805747 ,
135                                          0.0928785,
136                                          0.0625243 };
137
138   Double_t prior_pbpb[AliPID::kSPECIES] = {0.0609,
139                                            0.1064,
140                                            0.7152 ,
141                                            0.0442,
142                                            0.0733 };
143   cutPID->GetNegCut()->SetPriors(prior_pp);
144   cutPID->GetPosCut()->SetPriors(prior_pp);
145   cutPID->GetNegCut()->SetDetectors("ITS TPC TRD");
146   cutPID->GetPosCut()->SetDetectors("ITS TPC TRD");
147
148   if (readAOD) {
149     cutPID->GetNegCut()->SetAODmode(kTRUE);
150     cutPID->GetPosCut()->SetAODmode(kTRUE);
151   }
152   else {
153     cutPID->GetNegCut()->SetAODmode(kFALSE);
154     cutPID->GetPosCut()->SetAODmode(kFALSE);
155   }
156
157   cutPID->GetNegCut()->SetProbabilityCut(0);
158   cutPID->GetPosCut()->SetProbabilityCut(0);
159
160   switch(TMath::Abs(PDG)) {
161   case  310    : 
162     cutPID->GetNegCut()->SetParticleType(AliPID::kPion  , kTRUE);
163     cutPID->GetPosCut()->SetParticleType(AliPID::kPion  , kTRUE);
164     break;
165   case  3122   : 
166     cutPID->GetNegCut()->SetParticleType(AliPID::kPion  , kTRUE);
167     cutPID->GetPosCut()->SetParticleType(AliPID::kProton, kTRUE);
168     break;
169   case -3122   : 
170     cutPID->GetNegCut()->SetParticleType(AliPID::kProton, kTRUE);
171     cutPID->GetPosCut()->SetParticleType(AliPID::kPion  , kTRUE);
172     break;
173   default      : printf("UNDEFINED PID\n"); break;
174   }
175
176   Info("AliCFV0Task","CREATE MC KINE CUTS");
177   TObjArray* mcList = new TObjArray(0) ;
178   mcList->AddLast(mcKineCuts);
179   mcList->AddLast(mcGenCuts);
180
181   Info("AliCFV0Task","CREATE ACCEPTANCE CUTS");
182   TObjArray* accList = new TObjArray(0) ;
183   accList->AddLast(mcAccCuts);
184
185   Info("AliCFV0Task","CREATE RECONSTRUCTION CUTS");
186   TObjArray* recList = new TObjArray(0) ;
187   recList->AddLast(recKineCuts);
188   recList->AddLast(recQualityCuts);
189   recList->AddLast(recTopoCuts);
190
191   Info("AliCFV0Task","CREATE PID CUTS");
192   TObjArray* fPIDCutList = new TObjArray(0) ;
193   fPIDCutList->AddLast(cutPID);
194
195
196   //CREATE THE INTERFACE TO CORRECTION FRAMEWORK USED IN THE TASK
197   Info("AliCFV0Task","CREATE INTERFACE AND CUTS");
198   AliCFManager* man = new AliCFManager() ;
199   man->SetParticleContainer(container);
200   man->SetParticleCutsList (AliCFManager::kPartGenCuts,mcList);
201   man->SetParticleCutsList (AliCFManager::kPartAccCuts,accList);
202   man->SetParticleCutsList (AliCFManager::kPartRecCuts,recList);
203   man->SetParticleCutsList (AliCFManager::kPartSelCuts,fPIDCutList);
204
205   //CREATE THE TASK
206   Info("AliCFV0Task","CREATE TASK");
207   // create the task
208   AliCFV0Task *task = new AliCFV0Task("AliCFV0Task");
209   task->SetCFManager(man); //here is set the CF manager
210   if (!readAOD)      task->SetRebuildV0s(kTRUE);
211   task->SetV0PDG(PDG);
212
213   //SETUP THE ANALYSIS MANAGER TO READ INPUT CHAIN AND WRITE DESIRED OUTPUTS
214   Info("AliCFV0Task","CREATE ANALYSIS MANAGER");
215   // Make the analysis manager
216   AliAnalysisManager *mgr = new AliAnalysisManager("TestManager");
217
218   if (useGrid) mgr->SetAnalysisType(AliAnalysisManager::kGridAnalysis);
219   else mgr->SetAnalysisType(AliAnalysisManager::kLocalAnalysis);
220
221   AliMCEventHandler* mcHandler = new AliMCEventHandler();
222   mgr->SetMCtruthEventHandler(mcHandler);
223
224   AliInputEventHandler* dataHandler ;
225
226   if   (readAOD) dataHandler = new AliAODInputHandler();
227   else           dataHandler = new AliESDInputHandler();
228   mgr->SetInputEventHandler(dataHandler);
229
230   // Create and connect containers for input/output
231
232   //input data
233   AliAnalysisDataContainer *cinput0  = 
234     mgr->CreateContainer("cchain0",TChain::Class(),AliAnalysisManager::kInputContainer);
235   //slot 0 : default output tree (by default handled by AliAnalysisTaskSE)
236   AliAnalysisDataContainer *coutput0 = mgr->CreateContainer("ctree0", TTree::Class(),AliAnalysisManager::kOutputContainer,"outputV0.root");
237   // output histo (number of events processed)
238   AliAnalysisDataContainer *coutput1 = 
239     mgr->CreateContainer("chist0", TH1I::Class(),AliAnalysisManager::kOutputContainer,"outputV0.root");
240   // output Correction Framework Container (for acceptance & efficiency calculations)
241   AliAnalysisDataContainer *coutput2 = 
242     mgr->CreateContainer("ccontainer0", AliCFContainer::Class(),AliAnalysisManager::kOutputContainer,"outputV0.root");
243
244   cinput0->SetData(analysisChain);
245
246   mgr->AddTask(task);
247   mgr->ConnectInput (task,0,mgr->GetCommonInputContainer());
248   mgr->ConnectOutput(task,0,coutput0);
249   mgr->ConnectOutput(task,1,coutput1);
250   mgr->ConnectOutput(task,2,coutput2);
251
252  
253   Info("AliCFV0Task","READY TO RUN");
254   //RUN !!!
255   if (mgr->InitAnalysis()) {
256     mgr->PrintStatus();
257     mgr->StartAnalysis("local",analysisChain);
258   }
259
260   benchmark.Stop("AliCFV0Task");
261   benchmark.Show("AliCFV0Task");
262   
263   return kTRUE ;
264 }
265
266 void Load(Bool_t useGrid) {
267
268   //load the required aliroot libraries
269   gSystem->Load("libANALYSIS") ;
270   gSystem->Load("libANALYSISalice") ;
271   gSystem->Load("libCORRFW.so") ;
272
273   //compile online the task class
274   gSystem->SetIncludePath("-I. -I$ALICE_ROOT/include -I$ROOTSYS/include");
275   gROOT->LoadMacro("./AliCFV0Task.cxx+");
276 }