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