]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG3/vertexingHF/AliCFHeavyFlavourTask.C
New way to set input and output containers for AliAnalysisManager
[u/mrichter/AliRoot.git] / PWG3 / vertexingHF / AliCFHeavyFlavourTask.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 =  10.0 ;
6 const Int_t    mintrackrefsTPC = 2 ;
7 const Int_t    mintrackrefsITS = 3 ;
8 const Int_t    charge  = 1 ;
9 const Int_t    PDG = 421; 
10 const Int_t    minclustersTPC = 50 ;
11 //----------------------------------------------------
12
13 Bool_t AliCFHeavyFlavourTask()
14 {
15         
16         TBenchmark benchmark;
17         benchmark.Start("AliCFHeavyFlavourTask");
18         
19         AliLog::SetGlobalDebugLevel(0);
20         
21         //load the required libraries
22         Bool_t useParFiles=kFALSE;
23         gROOT->LoadMacro("$ALICE_ROOT/PWG3/vertexingHF/LoadLibraries.C");
24         LoadLibraries(useParFiles);
25
26         
27
28         TChain * analysisChain ;
29         TChain * analysisChainFriend ;
30         
31         //here put your input data path
32         printf("\n\nRunning on local file, please check the path\n\n");
33         analysisChain = new TChain("aodTree");
34         analysisChainFriend = new TChain("aodTree");
35         //for (Int_t i = 0; i<99; i++){
36         //TString fileName = "/Events/180002/";  // set here your own path!
37         //TString aodHFFileName = "/Events/180002/";  // set here your own path!
38         //fileName+=Form("%i/AliAOD.root",i);
39         //printf(Form("file = %s \n", fileName.Data()));
40         //aodHFFileName+=Form("%i/AliAOD.VertexingHF.root",i);
41         //analysisChain->Add(fileName);
42         //analysisChainFriend->Add(aodHFFileName);
43         analysisChain->Add("AliAOD.root");
44         analysisChainFriend->Add("AliAOD.VertexingHF.root");
45           //}
46                         
47         analysisChain->AddFriend(analysisChainFriend);
48         Info("AliCFHeavyFlavourTask",Form("CHAIN HAS %d ENTRIES",(Int_t)analysisChain->GetEntries()));
49         
50         //CONTAINER DEFINITION
51         Info("AliCFHeavyFlavourTask","SETUP CONTAINER");
52         //the sensitive variables (2 in this example), their indices
53         UInt_t ipt = 0;
54         UInt_t iy  = 1;
55         //Setting up the container grid... 
56         UInt_t nstep = 2 ; //number of selection steps MC 
57         const Int_t nvar   = 2 ; //number of variables on the grid:pt,y
58         const Int_t nbin1  = 8 ; //bins in pt
59         const Int_t nbin2  = 8 ; //bins in y 
60         
61         //arrays for the number of bins in each dimension
62         Int_t iBin[nvar];
63         iBin[0]=nbin1;
64         iBin[1]=nbin2;
65         
66         //arrays for lower bounds :
67         Double_t *binLim1=new Double_t[nbin1+1];
68         Double_t *binLim2=new Double_t[nbin2+1];
69         
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
74         //one "container" for MC
75         AliCFContainer* container = new AliCFContainer("container","container for tracks",nstep,nvar,iBin);
76
77         //setting the bin limits
78         container -> SetBinLimits(ipt,binLim1);
79         container -> SetBinLimits(iy,binLim2);
80         
81         //CREATE THE  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         //Particle-Level cuts:  
90         AliCFParticleGenCuts* mcGenCuts = new AliCFParticleGenCuts("mcGenCuts","MC particle generation cuts");
91         //mcGenCuts->SetRequireIsPrimary();
92         mcGenCuts->SetRequirePdgCode(PDG);
93         mcGenCuts->SetAODMC(1); //special flag for reading MC in AOD tree (important)
94         
95         // Rec-Level kinematic cuts
96         AliCFTrackKineCuts *recKineCuts = new AliCFTrackKineCuts("recKineCuts","rec-level kine cuts");
97         //   recKineCuts->SetPtRange(ptmin,ptmax);
98         //   recKineCuts->SetRapidityRange(ymin,ymax);
99         //   recKineCuts->SetChargeRec(charge);
100         
101         AliCFTrackQualityCuts *recQualityCuts = new AliCFTrackQualityCuts("recQualityCuts","rec-level quality cuts");
102         //recQualityCuts->SetStatus(AliESDtrack::kITSrefit);
103         
104         AliCFTrackIsPrimaryCuts *recIsPrimaryCuts = new AliCFTrackIsPrimaryCuts("recIsPrimaryCuts","rec-level isPrimary cuts");
105         //recIsPrimaryCuts->SetAODType(AliAODTrack::kPrimary);
106         
107         printf("CREATE MC KINE CUTS\n");
108         TObjArray* mcList = new TObjArray(0) ;
109         mcList->AddLast(mcKineCuts);
110         mcList->AddLast(mcGenCuts);
111         
112         printf("CREATE RECONSTRUCTION CUTS\n");
113         TObjArray* recList = new TObjArray(0) ;
114         recList->AddLast(recKineCuts);
115         recList->AddLast(recQualityCuts);
116         recList->AddLast(recIsPrimaryCuts);
117         
118         //CREATE THE INTERFACE TO CORRECTION FRAMEWORK USED IN THE TASK
119         printf("CREATE INTERFACE AND CUTS\n");
120         AliCFManager* man = new AliCFManager() ;
121         man->SetParticleContainer     (container);
122         man->SetParticleCutsList(0 , mcList); // MC
123         man->SetParticleCutsList(1 , recList); // AOD
124         
125         //CREATE THE TASK
126         printf("CREATE TASK\n");
127         // create the task
128         AliCFHeavyFlavourTask *task = new AliCFHeavyFlavourTask("AliCFHeavyFlavourTask");
129         task->SetCFManager(man); //here is set the CF manager
130         
131         //SETUP THE ANALYSIS MANAGER TO READ INPUT CHAIN AND WRITE DESIRED OUTPUTS
132         printf("CREATE ANALYSIS MANAGER\n");
133         // Make the analysis manager
134         AliAnalysisManager *mgr = new AliAnalysisManager("TestManager");
135         
136         mgr->SetAnalysisType(AliAnalysisManager::kLocalAnalysis);
137         
138         AliAODInputHandler* dataHandler = new AliAODInputHandler();
139         mgr->SetInputEventHandler(dataHandler);
140         
141         // Create and connect containers for input/output
142         
143         // ------ input data ------
144         AliAnalysisDataContainer *cinput0  = mgr->CreateContainer("cchain0",TChain::Class(),AliAnalysisManager::kInputContainer);
145         
146         // ----- output data -----
147         
148         //slot 0 : default output tree (by default handled by AliAnalysisTaskSE)
149         AliAnalysisDataContainer *coutput0 = mgr->CreateContainer("ctree0", TTree::Class(),AliAnalysisManager::kOutputContainer,"outputEta.root");
150         
151         //now comes user's output objects :
152         
153         // output TH1I for event counting
154         AliAnalysisDataContainer *coutput1 = mgr->CreateContainer("chist0", TH1I::Class(),AliAnalysisManager::kOutputContainer,"output.root");
155         // output Correction Framework Container (for acceptance & efficiency calculations)
156         AliAnalysisDataContainer *coutput2 = mgr->CreateContainer("ccontainer0", AliCFContainer::Class(),AliAnalysisManager::kOutputContainer,"output.root");
157         
158         cinput0->SetData(analysisChain);
159         
160         mgr->AddTask(task);
161         mgr->ConnectInput(task,0,mgr->GetCommonInputContainer());
162         mgr->ConnectOutput(task,0,coutput0);
163         mgr->ConnectOutput(task,1,coutput1);
164         mgr->ConnectOutput(task,2,coutput2);
165         
166         printf("READY TO RUN\n");
167         //RUN !!!
168         if (mgr->InitAnalysis()) {
169                 mgr->PrintStatus();
170                 mgr->StartAnalysis("local",analysisChain);
171         }
172         
173         benchmark.Stop("AliCFHeavyFlavourTask");
174         benchmark.Show("AliCFHeavyFlavourTask");
175         
176         return kTRUE ;
177 }
178