]> git.uio.no Git - u/mrichter/AliRoot.git/blob - CORRFW/test/AliCFRsnTask.C
53be9a2a2d8b1fe0fde62819fd9acaabf6c2dc4a
[u/mrichter/AliRoot.git] / CORRFW / test / AliCFRsnTask.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  = 0 ;
9 const Int_t      PDG = 313; 
10 const Int_t      minclustersTPC = 50 ;
11 const Double32_t nsigmavtx = 3. ;       //max track sigma to PVertex
12 //----------------------------------------------------
13
14 Bool_t AliCFRsnTask(
15                     const Bool_t useGrid = 1,
16                     const Bool_t readAOD = 0,
17                     const char * kTagXMLFile="wn.xml", // XML file containing tags
18                     )
19 {
20   
21   TBenchmark benchmark;
22   benchmark.Start("AliRsnTask");
23
24   AliLog::SetGlobalDebugLevel(0);
25
26   Load(useGrid) ; //load the required libraries
27
28   TChain * analysisChain ;
29
30   if (useGrid) { //data located on AliEn
31     TGrid::Connect("alien://") ;    //  Create an AliRunTagCuts and an AliEventTagCuts Object and impose some selection criteria
32     AliRunTagCuts      *runCuts   = new AliRunTagCuts(); 
33     AliEventTagCuts    *eventCuts = new AliEventTagCuts(); 
34     AliLHCTagCuts      *lhcCuts   = new AliLHCTagCuts(); 
35     AliDetectorTagCuts *detCuts   = new AliDetectorTagCuts(); 
36     eventCuts->SetMultiplicityRange(0,20000);
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     //  Create a new esd chain and assign the chain that is returned by querying the tags
45     analysisChain = tagAna->QueryTags(runCuts,lhcCuts,detCuts,eventCuts); 
46   }
47   else {// local data
48     //here put your input data path
49     if (readAOD) {
50       analysisChain = new TChain("aodTree");
51       //analysisChain->Add("AliAOD.root");
52       analysisChain->Add("/home/vernet/Data/LHC08b2/300000/001/AliAOD.root");
53       analysisChain->Add("/home/vernet/Data/LHC08b2/300000/002/AliAOD.root");
54       analysisChain->Add("/home/vernet/Data/LHC08b2/300000/003/AliAOD.root");
55       analysisChain->Add("/home/vernet/Data/LHC08b2/300000/004/AliAOD.root");
56     }
57     else {
58       analysisChain = new TChain("esdTree");
59       //analysisChain->Add("AliESDs.root");
60       analysisChain->Add("/home/vernet/Data/LHC08b2/300000/001/AliESDs.root");
61       analysisChain->Add("/home/vernet/Data/LHC08b2/300000/002/AliESDs.root");
62       analysisChain->Add("/home/vernet/Data/LHC08b2/300000/003/AliESDs.root");
63       analysisChain->Add("/home/vernet/Data/LHC08b2/300000/004/AliESDs.root");
64     }
65   }
66   
67   
68   Info("AliCFRsnTask",Form("CHAIN HAS %d ENTRIES",(Int_t)analysisChain->GetEntries()));
69
70   //CONTAINER DEFINITION
71   Info("AliCFRsnTask","SETUP CONTAINER");
72   //the sensitive variables (2 in this example), their indices
73   Int_t ipt = 0;
74   Int_t iy  = 1;
75   //Setting up the container grid... 
76   Int_t nstep = 4 ; //number of selection steps MC 
77   const Int_t nvar   = 2 ; //number of variables on the grid:pt,y
78   const Int_t nbin1  = 8 ; //bins in pt
79   const Int_t nbin2  = 8 ; //bins in y 
80   //arrays for the number of bins in each dimension
81   const Int_t iBin[nvar] ={nbin1,nbin2};
82   //arrays for lower bounds :
83   Double_t binLim1[nbin1+1];
84   Double_t binLim2[nbin2+1];
85   //values for bin lower bounds
86   for(Int_t i=0; i<=nbin1; i++) binLim1[i]=(Double_t)ptmin + (ptmax-ptmin)/nbin1*(Double_t)i ; 
87   for(Int_t i=0; i<=nbin2; i++) binLim2[i]=(Double_t)ymin  + (ymax-ymin)  /nbin2*(Double_t)i ;
88   //one "container" for MC
89   AliCFContainer* container = new AliCFContainer("container","container for tracks",nstep,nvar,iBin);
90   //setting the bin limits
91   container -> SetBinLimits(ipt,binLim1);
92   container -> SetBinLimits(iy,binLim2);
93
94
95   //CREATE THE  CUTS -----------------------------------------------
96   //Particle-Level cuts:  
97
98   // Gen-Level kinematic cuts
99   AliCFTrackKineCuts *mcKineCuts = new AliCFTrackKineCuts("mcKineCuts","MC-level kinematic cuts");
100   mcKineCuts->SetPtRange(ptmin,ptmax);
101   mcKineCuts->SetRapidityRange(ymin,ymax);
102   mcKineCuts->SetChargeMC(charge);
103
104   AliCFParticleGenCuts* mcGenCuts = new AliCFParticleGenCuts("mcGenCuts","MC particle generation cuts");
105   //  mcGenCuts->SetRequireIsPrimary();
106   mcGenCuts->SetRequirePdgCode(PDG);
107
108   //Acceptance Cuts
109   AliCFPairAcceptanceCuts *mcAccCuts = new AliCFPairAcceptanceCuts("mcAccCuts","MC acceptance cuts");
110   mcAccCuts->SetMinNHitITS(mintrackrefsITS,mintrackrefsITS);
111   mcAccCuts->SetMinNHitTPC(mintrackrefsTPC,mintrackrefsTPC);
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
119   AliCFPairQualityCuts *recQualityCuts = new AliCFPairQualityCuts("recQualityCuts","rec-level quality cuts");
120   if (!readAOD) recQualityCuts->SetMinNClusterTPC(minclustersTPC,minclustersTPC);
121   recQualityCuts->SetStatus(AliESDtrack::kTPCrefit & AliESDtrack::kITSrefit, 
122                             AliESDtrack::kTPCrefit & AliESDtrack::kITSrefit) ;
123
124
125
126   AliCFPairIsPrimaryCuts *recIsPrimaryCuts = new AliCFPairIsPrimaryCuts("recIsPrimaryCuts","rec-level isPrimary cuts");
127   if (readAOD) recIsPrimaryCuts->SetAODType(AliAODTrack::kPrimary,AliAODTrack::kPrimary);
128   else         recIsPrimaryCuts->SetMaxNSigmaToVertex(nsigmavtx,nsigmavtx);
129
130   AliCFPairPidCut* cutPID = new AliCFPairPidCut("cutPID","ESD_PID") ;
131   Double_t prior_pp[AliPID::kSPECIES] = {0.0244519,
132                                          0.0143988,
133                                          0.805747 ,
134                                          0.0928785,
135                                          0.0625243 };
136
137   Double_t prior_pbpb[AliPID::kSPECIES] = {0.0609,
138                                            0.1064,
139                                            0.7152 ,
140                                            0.0442,
141                                            0.0733 };
142
143
144   cutPID->SetPriors(prior_pbpb);
145   cutPID->SetProbabilityCut(0.,0.);
146   cutPID->SetDetectors("TPC ITS TOF TRD","TPC ITS TOF TRD");
147   if (readAOD) cutPID->SetAODmode(kTRUE);
148   else         cutPID->SetAODmode(kFALSE);
149   
150   switch(PDG) {
151   case  -313  : cutPID->SetParticleType(AliPID::kKaon  ,kTRUE,AliPID::kPion  ,kTRUE); break;
152   case   313  : cutPID->SetParticleType(AliPID::kPion  ,kTRUE,AliPID::kKaon  ,kTRUE); break;
153   case   333  : cutPID->SetParticleType(AliPID::kKaon  ,kTRUE,AliPID::kKaon  ,kTRUE); break;
154   case  3124  : cutPID->SetParticleType(AliPID::kKaon  ,kTRUE,AliPID::kProton,kTRUE); break;
155   case -3124  : cutPID->SetParticleType(AliPID::kProton,kTRUE,AliPID::kKaon  ,kTRUE); break;
156   default     : printf("UNDEFINED PID\n"); break;
157   }
158   //cutPID->SetQAOn(kTRUE);
159
160   Info("AliCFRsnTask","CREATE MC KINE CUTS");
161   TObjArray* mcList = new TObjArray(0) ;
162   mcList->AddLast(mcKineCuts);
163   mcList->AddLast(mcGenCuts);
164
165   Info("AliCFRsnTask","CREATE ACCEPTANCE CUTS");
166   TObjArray* accList = new TObjArray(0) ;
167   accList->AddLast(mcAccCuts);
168
169   Info("AliCFRsnTask","CREATE RECONSTRUCTION CUTS");
170   TObjArray* recList = new TObjArray(0) ;
171   recList->AddLast(recKineCuts);
172   recList->AddLast(recQualityCuts);
173   recList->AddLast(recIsPrimaryCuts);
174
175   Info("AliCFRsnTask","CREATE PID CUTS");
176   TObjArray* fPIDCutList = new TObjArray(0) ;
177   fPIDCutList->AddLast(cutPID);
178
179
180   //CREATE THE INTERFACE TO CORRECTION FRAMEWORK USED IN THE TASK
181   Info("AliCFRsnTask","CREATE INTERFACE AND CUTS");
182   AliCFManager* man = new AliCFManager() ;
183   man->SetParticleContainer     (container);
184   man->SetParticleCutsList(AliCFManager::kPartGenCuts,mcList);
185   man->SetParticleCutsList(AliCFManager::kPartAccCuts,accList);
186   man->SetParticleCutsList(AliCFManager::kPartRecCuts,recList);
187   man->SetParticleCutsList(AliCFManager::kPartSelCuts,fPIDCutList);
188
189
190   //CREATE THE TASK
191   Info("AliCFRsnTask","CREATE TASK");
192   // create the task
193   AliCFRsnTask *task = new AliCFRsnTask("AliRsnTask");
194   task->SetCFManager(man); //here is set the CF manager
195   task->SetRsnPDG(PDG);
196
197
198   //SETUP THE ANALYSIS MANAGER TO READ INPUT CHAIN AND WRITE DESIRED OUTPUTS
199   Info("AliCFRsnTask","CREATE ANALYSIS MANAGER");
200   // Make the analysis manager
201   AliAnalysisManager *mgr = new AliAnalysisManager("TestManager");
202
203   if (useGrid) mgr->SetAnalysisType(AliAnalysisManager::kGridAnalysis);
204   else mgr->SetAnalysisType(AliAnalysisManager::kLocalAnalysis);
205
206   AliMCEventHandler* mcHandler = new AliMCEventHandler();
207   mgr->SetMCtruthEventHandler(mcHandler);
208
209   AliInputEventHandler* dataHandler ;
210   if (readAOD) dataHandler = new AliAODInputHandler();
211   else         dataHandler = new AliESDInputHandler();
212   mgr->SetInputEventHandler(dataHandler);
213   
214
215
216   // Create and connect containers for input/output
217
218   //input data
219   AliAnalysisDataContainer *cinput0  = mgr->CreateContainer("cchain0",TChain::Class(),AliAnalysisManager::kInputContainer);
220
221   //slot 0 : default output tree (by default handled by AliAnalysisTaskSE)
222   AliAnalysisDataContainer *coutput0 = mgr->CreateContainer("ctree0", TTree::Class(),AliAnalysisManager::kOutputContainer,"output_rsn.root");
223
224   // output histo (number of events processed)
225   AliAnalysisDataContainer *coutput1 = mgr->CreateContainer("chist0", TH1I::Class(),AliAnalysisManager::kOutputContainer,"output_rsn.root");
226   // output Correction Framework Container (for acceptance & efficiency calculations)
227   AliAnalysisDataContainer *coutput2 = mgr->CreateContainer("ccontainer0", AliCFContainer::Class(),AliAnalysisManager::kOutputContainer,"output_rsn.root");
228
229   cinput0->SetData(analysisChain);
230
231   mgr->AddTask(task);
232   mgr->ConnectInput (task,0,mgr->GetCommonInputContainer());
233   mgr->ConnectOutput(task,0,coutput0);
234   mgr->ConnectOutput(task,1,coutput1);
235   mgr->ConnectOutput(task,2,coutput2);
236  
237
238   Info("AliCFRsnTask","READY TO RUN");
239   //RUN !!!
240   if (mgr->InitAnalysis()) {
241     mgr->PrintStatus();
242     mgr->StartAnalysis("local",analysisChain);
243   }
244
245   benchmark.Stop("AliRsnTask");
246   benchmark.Show("AliRsnTask");
247   
248   return kTRUE ;
249 }
250
251 void Load(Bool_t useGrid) {
252   //load the required aliroot libraries
253   gSystem->Load("libANALYSIS") ;
254   gSystem->Load("libANALYSISalice") ;
255   gSystem->Load("libPWG2resonances");
256   gSystem->Load("libCORRFW");
257
258   gSystem->AddIncludePath("-I$ALICE_ROOT/PWG2/RESONANCES");
259   gROOT->LoadMacro("./AliCFRsnTask.cxx+");
260 }