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