]> git.uio.no Git - u/mrichter/AliRoot.git/blob - CORRFW/test/AliCFSingleTrackTask.C
add documentation and example macro
[u/mrichter/AliRoot.git] / CORRFW / test / AliCFSingleTrackTask.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  = 1 ;
9 const Int_t    PDG = 2212; 
10 const Int_t    minclustersTPC = 50 ;
11 //----------------------------------------------------
12
13 Bool_t AliCFSingleTrackTask(
14                             const Bool_t useGrid = 1,
15                             const char * kTagXMLFile="wn.xml", // XML file containing tags
16                             Long64_t     nentries=TChain::kBigNumber
17                             )
18 {
19   
20   TBenchmark benchmark;
21   benchmark.Start("AliSingleTrackTask");
22
23   AliLog::SetGlobalDebugLevel(0);
24
25   Load() ; //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("AliCFSingleTrackTask",Form("CHAIN HAS %d ENTRIES",(Int_t)analysisChain->GetEntries()));
54
55   //CONTAINER DEFINITION
56   Info("AliCFSingleTrackTask","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   AliCFAcceptanceCuts *mcAccCuts = new AliCFAcceptanceCuts("mcAccCuts","MC acceptance cuts");
95   mcAccCuts->SetMinNHitITS(mintrackrefsITS);
96   mcAccCuts->SetMinNHitTPC(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   // QA histograms for rec-level kinematic cuts
104   recKineCuts->SetQAOn(kTRUE);
105
106   AliCFTrackQualityCuts *recQualityCuts = new AliCFTrackQualityCuts("recQualityCuts","rec-level quality cuts");
107   recQualityCuts->SetMinNClusterTPC(minclustersTPC);
108   recQualityCuts->SetRequireITSRefit(kTRUE);
109   // QA histograms for rec-level quality cuts
110   recQualityCuts->SetQAOn(kTRUE);
111
112   AliCFTrackIsPrimaryCuts *recIsPrimaryCuts = new AliCFTrackIsPrimaryCuts("recIsPrimaryCuts","rec-level isPrimary cuts");
113   recIsPrimaryCuts->SetMaxNSigmaToVertex(3);
114   // QA histograms for rec-level primary-check cuts
115   recIsPrimaryCuts->SetQAOn(kTRUE);
116
117   AliCFTrackCutPid* cutPID = new AliCFTrackCutPid("cutPID","ESD_PID") ;
118   Double_t prior[AliPID::kSPECIES] = {0.0244519,
119                                       0.0143988,
120                                       0.805747 ,
121                                       0.0928785,
122                                       0.0625243 };
123   cutPID->SetPriors(prior);
124   cutPID->SetProbabilityCut(0.0);
125   cutPID->SetDetectors("TPC TOF");
126   switch(TMath::Abs(PDG)) {
127   case 11   : cutPID->SetParticleType(AliPID::kElectron, kTRUE); break;
128   case 13   : cutPID->SetParticleType(AliPID::kMuon    , kTRUE); break;
129   case 211  : cutPID->SetParticleType(AliPID::kPion    , kTRUE); break;
130   case 321  : cutPID->SetParticleType(AliPID::kKaon    , kTRUE); break;
131   case 2212 : cutPID->SetParticleType(AliPID::kProton  , kTRUE); break;
132   default   : printf("UNDEFINED PID\n"); break;
133   }
134   cutPID->SetQAOn(kTRUE);
135
136   printf("CREATE MC KINE CUTS\n");
137   TObjArray* mcList = new TObjArray(0) ;
138   mcList->AddLast(mcKineCuts);
139   mcList->AddLast(mcGenCuts);
140
141   printf("CREATE ACCEPTANCE CUTS\n");
142   TObjArray* accList = new TObjArray(0) ;
143   accList->AddLast(mcAccCuts);
144
145   printf("CREATE RECONSTRUCTION CUTS\n");
146   TObjArray* recList = new TObjArray(0) ;
147   recList->AddLast(recKineCuts);
148   recList->AddLast(recQualityCuts);
149   recList->AddLast(recIsPrimaryCuts);
150
151   printf("CREATE PID CUTS\n");
152   TObjArray* fPIDCutList = new TObjArray(0) ;
153   fPIDCutList->AddLast(cutPID);
154
155
156   //CREATE THE INTERFACE TO CORRECTION FRAMEWORK USED IN THE TASK
157   printf("CREATE INTERFACE AND CUTS\n");
158   AliCFManager* man = new AliCFManager() ;
159   man->SetParticleContainer     (container);
160   man->SetParticleCutsList(AliCFManager::kPartGenCuts,mcList);
161   man->SetParticleCutsList(AliCFManager::kPartAccCuts,accList);
162   man->SetParticleCutsList(AliCFManager::kPartRecCuts,recList);
163   man->SetParticleCutsList(AliCFManager::kPartSelCuts,fPIDCutList);
164
165
166   //CREATE THE TASK
167   printf("CREATE TASK\n");
168   // create the task
169   AliCFSingleTrackTask *task = new AliCFSingleTrackTask("AliSingleTrackTask");
170   task->SetCFManager(man); //here is set the CF manager
171
172
173   //SETUP THE ANALYSIS MANAGER TO READ INPUT CHAIN AND WRITE DESIRED OUTPUTS
174   printf("CREATE ANALYSIS MANAGER\n");
175   // Make the analysis manager
176   AliAnalysisManager *mgr = new AliAnalysisManager("TestManager");
177
178   // Create and connect containers for input/output
179
180   //input data
181   AliAnalysisDataContainer *cinput0  = mgr->CreateContainer("cchain0",TChain::Class(),AliAnalysisManager::kInputContainer);
182   // output histo (number of events processed)
183   AliAnalysisDataContainer *coutput0 = mgr->CreateContainer("chist0", TH1I::Class(),AliAnalysisManager::kOutputContainer,"output.root");
184   // output Correction Framework Container (for acceptance & efficiency calculations)
185   AliAnalysisDataContainer *coutput1 = mgr->CreateContainer("ccontainer0", AliCFContainer::Class(),AliAnalysisManager::kOutputContainer,"output.root");
186   // output QA histograms 
187   AliAnalysisDataContainer *coutput2 = mgr->CreateContainer("clist0", TList::Class(),AliAnalysisManager::kOutputContainer,"output.root");
188
189   mgr->AddTask(task);
190   mgr->ConnectInput(task,0,cinput0);
191   mgr->ConnectOutput(task,0,coutput0);
192   mgr->ConnectOutput(task,1,coutput1);
193   mgr->ConnectOutput(task,2,coutput2);
194   cinput0->SetData(analysisChain);
195  
196   //NEW INTERFACE TO MC INFORMATION
197   AliMCEventHandler* mcHandler = new AliMCEventHandler();
198   mgr->SetMCtruthEventHandler(mcHandler);
199
200   printf("READY TO RUN\n");
201   //RUN !!!
202   if (mgr->InitAnalysis()) {
203     mgr->PrintStatus();
204     mgr->StartAnalysis("local",analysisChain);
205   }
206
207   benchmark.Stop("AliSingleTrackTask");
208   benchmark.Show("AliSingleTrackTask");
209   
210   return kTRUE ;
211 }
212
213 void Load() {
214   //remove this file which can cause problems
215   gSystem->Exec("rm $ALICE_ROOT/ANALYSIS/AliAnalysisSelector_cxx.so");
216
217   //load the required aliroot libraries
218   gSystem->Load("libANALYSIS") ;
219   gSystem->Load("libANALYSISRL") ;
220   gSystem->Load("$ALICE_ROOT/CORRFW/libCORRFW.so") ;
221
222   //compile online the task class
223   gSystem->SetIncludePath("-I. -I$ALICE_ROOT/include -I$ROOTSYS/include");
224   gROOT->LoadMacro("./AliCFSingleTrackTask.cxx+");
225 }