]> git.uio.no Git - u/mrichter/AliRoot.git/blob - CORRFW/test/AliCFRsnTask.C
e211eb5732c32e65df67841b4ed6cd7df1886c28
[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 //   AliCFPairKineCuts *recKineCuts = new AliCFPairKineCuts("recKineCuts","rec-level kine cuts");
101   recKineCuts->SetPtRange(ptmin,ptmax);
102   recKineCuts->SetRapidityRange(ymin,ymax);
103   recKineCuts->SetChargeRec(charge);
104
105   AliCFPairQualityCuts *recQualityCuts = new AliCFPairQualityCuts("recQualityCuts","rec-level quality cuts");
106   recQualityCuts->SetMinNClusterTPC(minclustersTPC,minclustersTPC);
107   recQualityCuts->SetRequireITSRefit(kTRUE,kTRUE);
108
109   AliCFPairIsPrimaryCuts *recIsPrimaryCuts = new AliCFPairIsPrimaryCuts("recIsPrimaryCuts","rec-level isPrimary cuts");
110   recIsPrimaryCuts->SetMaxNSigmaToVertex(nsigmavtx,nsigmavtx);
111
112   AliCFPairPidCut* cutPID = new AliCFPairPidCut("cutPID","ESD_PID") ;
113   Double_t prior_pp[AliPID::kSPECIES] = {0.0244519,
114                                          0.0143988,
115                                          0.805747 ,
116                                          0.0928785,
117                                          0.0625243 };
118
119   Double_t prior_pbpb[AliPID::kSPECIES] = {0.0609,
120                                            0.1064,
121                                            0.7152 ,
122                                            0.0442,
123                                            0.0733 };
124
125
126   cutPID->SetPriors(prior_pbpb);
127   cutPID->SetProbabilityCut(0.,0.);
128   cutPID->SetDetectors("TPC ITS TOF TRD","TPC ITS TOF TRD");
129   switch(PDG) {
130   case  -313  : cutPID->SetParticleType(AliPID::kKaon,kTRUE,AliPID::kPion,kTRUE); break;
131   case   313  : cutPID->SetParticleType(AliPID::kPion,kTRUE,AliPID::kKaon,kTRUE); break;
132   case   333  : cutPID->SetParticleType(AliPID::kKaon,kTRUE,AliPID::kKaon,kTRUE); break;
133   case  3124  : cutPID->SetParticleType(AliPID::kPion,kTRUE,AliPID::kKaon,kTRUE); break;
134   case -3124  : cutPID->SetParticleType(AliPID::kPion,kTRUE,AliPID::kKaon,kTRUE); break;
135   default     : printf("UNDEFINED PID\n"); break;
136   }
137 //   cutPID->SetQAOn(kTRUE);
138
139   Info("AliCFRsnTask","CREATE MC KINE CUTS");
140   TObjArray* mcList = new TObjArray(0) ;
141   mcList->AddLast(mcKineCuts);
142   mcList->AddLast(mcGenCuts);
143
144   Info("AliCFRsnTask","CREATE ACCEPTANCE CUTS");
145   TObjArray* accList = new TObjArray(0) ;
146   accList->AddLast(mcAccCuts);
147
148   Info("AliCFRsnTask","CREATE RECONSTRUCTION CUTS");
149   TObjArray* recList = new TObjArray(0) ;
150   recList->AddLast(recKineCuts);
151   recList->AddLast(recQualityCuts);
152   recList->AddLast(recIsPrimaryCuts);
153
154   Info("AliCFRsnTask","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("AliCFRsnTask","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
169   //CREATE THE TASK
170   Info("AliCFRsnTask","CREATE TASK");
171   // create the task
172   AliCFRsnTask *task = new AliCFRsnTask("AliRsnTask");
173   task->SetCFManager(man); //here is set the CF manager
174
175
176   //SETUP THE ANALYSIS MANAGER TO READ INPUT CHAIN AND WRITE DESIRED OUTPUTS
177   Info("AliCFRsnTask","CREATE ANALYSIS MANAGER");
178   // Make the analysis manager
179   AliAnalysisManager *mgr = new AliAnalysisManager("TestManager");
180
181   // Create and connect containers for input/output
182
183   //input data
184   AliAnalysisDataContainer *cinput0  = mgr->CreateContainer("cchain0",TChain::Class(),AliAnalysisManager::kInputContainer);
185   // output histo (number of events processed)
186   AliAnalysisDataContainer *coutput0 = mgr->CreateContainer("chist0", TH1I::Class(),AliAnalysisManager::kOutputContainer,"output_rsn.root");
187   // output Correction Framework Container (for acceptance & efficiency calculations)
188   AliAnalysisDataContainer *coutput1 = mgr->CreateContainer("ccontainer0", AliCFContainer::Class(),AliAnalysisManager::kOutputContainer,"output_rsn.root");
189   // output QA histograms 
190
191   mgr->AddTask(task);
192   mgr->ConnectInput(task,0,cinput0);
193   mgr->ConnectOutput(task,0,coutput0);
194   mgr->ConnectOutput(task,1,coutput1);
195   cinput0->SetData(analysisChain);
196  
197   //NEW INTERFACE TO MC INFORMATION
198   AliMCEventHandler* mcHandler = new AliMCEventHandler();
199   mgr->SetMCtruthEventHandler(mcHandler);
200
201   Info("AliCFRsnTask","READY TO RUN");
202   //RUN !!!
203   if (mgr->InitAnalysis()) {
204     mgr->PrintStatus();
205     mgr->StartAnalysis("local",analysisChain);
206   }
207
208   benchmark.Stop("AliRsnTask");
209   benchmark.Show("AliRsnTask");
210   
211   return kTRUE ;
212 }
213
214 void Load(Bool_t useGrid) {
215   //remove this file which can cause problems
216   if (!useGrid) 
217     gSystem->Exec("rm $ALICE_ROOT/ANALYSIS/AliAnalysisSelector_cxx.so");
218   
219   //load the required aliroot libraries
220   gSystem->Load("libANALYSIS") ;
221   
222   gSystem->SetIncludePath("-I. -I$ALICE_ROOT/include -I$ROOTSYS/include -I$ALICE_ROOT/PWG2/RESONANCES -I$ALICE_ROOT/CORRFW");
223
224   if (!useGrid) {
225     //load local CF library
226     gSystem->Load("libCORRFW");
227     //load PWG2 library
228     gSystem->Load("libPWG2resonances");
229   }
230   else {
231     //compile online the task class
232     gSystem->Exec("cp $ALICE_ROOT/PWG2/RESONANCES/AliRsnDaughter.cxx .");
233     gROOT->LoadMacro("./AliRsnDaughter.cxx+");
234     gSystem->Exec("tar zxvf cf_classes.tgz");
235     gROOT->LoadMacro("./AliCFFrame.cxx+");
236     gROOT->LoadMacro("./AliCFVGrid.cxx+");
237     gROOT->LoadMacro("./AliCFGrid.cxx+");
238     gROOT->LoadMacro("./AliCFGridSparse.cxx+");
239     gROOT->LoadMacro("./AliCFContainer.cxx+");
240     gROOT->LoadMacro("./AliCFCutBase.cxx+");
241     gROOT->LoadMacro("./AliCFTrackKineCuts.cxx+");
242     gROOT->LoadMacro("./AliCFTrackIsPrimaryCuts.cxx+");
243     gROOT->LoadMacro("./AliCFTrackQualityCuts.cxx+");
244     gROOT->LoadMacro("./AliCFParticleGenCuts.cxx+");
245     gROOT->LoadMacro("./AliCFAcceptanceCuts.cxx+");
246     gROOT->LoadMacro("./AliCFTrackCutPid.cxx+");
247     gROOT->LoadMacro("./AliCFManager.cxx+");
248     gSystem->Exec("tar zxvf cf_rsn.tgz");
249     gROOT->LoadMacro("./AliCFPair.cxx+");
250     gROOT->LoadMacro("./AliCFPairAcceptanceCuts.cxx+");
251     gROOT->LoadMacro("./AliCFPairQualityCuts.cxx+");
252     gROOT->LoadMacro("./AliCFPairIsPrimaryCuts.cxx+");
253     gROOT->LoadMacro("./AliCFPairPidCut.cxx+");
254   }
255   gROOT->LoadMacro("./AliCFRsnTask.cxx+");
256 }