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