]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/AliAnalysisTaskProtonsQA.cxx
new methods and possiility to run all methods in same run
[u/mrichter/AliRoot.git] / PWG2 / AliAnalysisTaskProtonsQA.cxx
1 #include "TChain.h"
2 #include "TTree.h"
3 #include "TString.h"
4 #include "TList.h"
5 #include "TH2F.h"
6 #include "TH1I.h"
7 #include "TF1.h"
8 #include "TCanvas.h"
9 #include "TLorentzVector.h"
10
11 #include "AliAnalysisTask.h"
12 #include "AliAnalysisManager.h"
13
14 #include "AliESDEvent.h"
15 #include "AliESDInputHandler.h"
16 #include "AliMCEventHandler.h"
17 #include "AliMCEvent.h"
18 #include "AliStack.h"
19
20 #include "PWG2spectra/SPECTRA/AliProtonQAAnalysis.h"
21 #include "AliAnalysisTaskProtonsQA.h"
22
23 // Analysis task used for the QA of the (anti)proton analysis
24 // Author: Panos Cristakoglou
25
26 ClassImp(AliAnalysisTaskProtonsQA)
27
28 //________________________________________________________________________ 
29 AliAnalysisTaskProtonsQA::AliAnalysisTaskProtonsQA()
30   : AliAnalysisTask(), fESD(0), fMC(0),
31     fList0(0), fList1(0), fList2(0), fList3(0), fList4(0), fList5(0),
32     fAnalysis(0) {
33     //Dummy constructor
34 }
35
36 //________________________________________________________________________
37 AliAnalysisTaskProtonsQA::AliAnalysisTaskProtonsQA(const char *name) 
38 : AliAnalysisTask(name, ""), fESD(0), fMC(0),
39   fList0(0), fList1(0), fList2(0), fList3(0), fList4(0), fList5(0),
40   fAnalysis(0) {
41   // Constructor
42
43   // Define input and output slots here
44   // Input slot #0 works with a TChain
45   DefineInput(0, TChain::Class());
46   // Output slot #0 writes into a TList container
47   DefineOutput(0, TList::Class());
48   DefineOutput(1, TList::Class());
49   DefineOutput(2, TList::Class());
50   DefineOutput(3, TList::Class());
51   DefineOutput(4, TList::Class());
52   DefineOutput(5, TList::Class());
53 }
54
55 //________________________________________________________________________
56 void AliAnalysisTaskProtonsQA::ConnectInputData(Option_t *) {
57   // Connect ESD or AOD here
58   // Called once
59
60   TTree* tree = dynamic_cast<TTree*> (GetInputData(0));
61   if (!tree) {
62     Printf("ERROR: Could not read chain from input slot 0");
63   } else {
64     AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
65      
66     if (!esdH) {
67       Printf("ERROR: Could not get ESDInputHandler");
68     } else
69       fESD = esdH->GetEvent();
70   }
71
72   AliMCEventHandler* mcH = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
73   if (!mcH) {
74     Printf("ERROR: Could not retrieve MC event handler");
75   }
76   else
77     fMC = mcH->MCEvent();
78 }
79
80 //________________________________________________________________________
81 void AliAnalysisTaskProtonsQA::CreateOutputObjects() {
82   // Create histograms
83   // Called once
84   //Prior probabilities
85   Double_t partFrac[5] = {0.01, 0.01, 0.85, 0.10, 0.05};
86   
87   //proton analysis object
88   fAnalysis = new AliProtonQAAnalysis();
89   fAnalysis->SetRunMCAnalysis();
90   fAnalysis->SetRunEfficiencyAnalysis();
91   //fAnalysis->SetMCProcessId(4);//4: weak decay - 13: hadronic interaction
92   //fAnalysis->SetMotherParticlePDGCode(3122);//3122: Lambda
93
94   //Use of TPConly tracks
95   fAnalysis->SetQAYPtBins(10, -0.5, 0.5, 12, 0.5, 0.9); //TPC only
96   fAnalysis->UseTPCOnly();
97   fAnalysis->SetMinTPCClusters(50);
98   fAnalysis->SetMaxChi2PerTPCCluster(2.5);
99   fAnalysis->SetMaxCov11(2.0);
100   fAnalysis->SetMaxCov22(2.0);
101   fAnalysis->SetMaxCov33(0.5);
102   fAnalysis->SetMaxCov44(0.5);
103   fAnalysis->SetMaxCov55(2.0);
104   fAnalysis->SetMaxSigmaToVertexTPC(2.0);
105   //fAnalysis->SetMaxDCAXYTPC(2.0);
106   //fAnalysis->SetMaxDCAZTPC(2.0);
107   fAnalysis->SetTPCpid();
108
109   //Use of HybridTPC tracks
110   /*fAnalysis->SetQAYPtBins(10, -0.5, 0.5, 12, 0.5, 0.9); //HybridTPC
111   fAnalysis->UseHybridTPC();
112   fAnalysis->SetMinTPCClusters(50);
113   fAnalysis->SetMaxChi2PerTPCCluster(2.5);
114   fAnalysis->SetMaxCov11(2.0);
115   fAnalysis->SetMaxCov22(2.0);
116   fAnalysis->SetMaxCov33(0.5);
117   fAnalysis->SetMaxCov44(0.5);
118   fAnalysis->SetMaxCov55(2.0);
119   fAnalysis->SetMaxSigmaToVertexTPC(2.0);
120   //fAnalysis->SetMaxDCAXY(2.0);
121   //fAnalysis->SetMaxDCAZ(2.0);
122   fAnalysis->SetTPCpid();
123   fAnalysis->SetPointOnITSLayer1();
124   fAnalysis->SetPointOnITSLayer2();
125   fAnalysis->SetPointOnITSLayer3();
126   fAnalysis->SetPointOnITSLayer4();
127   fAnalysis->SetPointOnITSLayer5();
128   fAnalysis->SetPointOnITSLayer6();
129   fAnalysis->SetMinITSClusters(5);*/
130
131   //Combined tracking
132   /*fAnalysis->SetQAYPtBins(20, -1.0, 1.0, 27, 0.4, 3.1); //combined tracking
133   fAnalysis->SetMinTPCClusters(50);
134   fAnalysis->SetMaxChi2PerTPCCluster(3.5);
135   fAnalysis->SetMaxCov11(2.0);
136   fAnalysis->SetMaxCov22(2.0);
137   fAnalysis->SetMaxCov33(0.5);
138   fAnalysis->SetMaxCov44(0.5);
139   fAnalysis->SetMaxCov55(2.0);
140   fAnalysis->SetMaxSigmaToVertex(2.0);
141   //fAnalysis->SetMaxDCAXY(2.0);
142   //fAnalysis->SetMaxDCAZ(2.0);
143   fAnalysis->SetTPCRefit();
144   fAnalysis->SetPointOnITSLayer1();
145   fAnalysis->SetPointOnITSLayer2();
146   fAnalysis->SetPointOnITSLayer3();
147   fAnalysis->SetPointOnITSLayer4();
148   fAnalysis->SetPointOnITSLayer5();
149   fAnalysis->SetPointOnITSLayer6();
150   fAnalysis->SetMinITSClusters(1);
151   fAnalysis->SetITSRefit();
152   fAnalysis->SetESDpid();*/
153
154   fAnalysis->SetPriorProbabilities(partFrac);
155
156   fList0 = new TList();
157   fList0 = fAnalysis->GetGlobalQAList();
158
159   fList1 = new TList();
160   fList1 = fAnalysis->GetPDGList();
161
162   fList2 = new TList();
163   fList2 = fAnalysis->GetMCProcessesList();
164
165   fList3 = new TList();
166   fList3 = fAnalysis->GetAcceptedCutList();
167
168   fList4 = new TList();
169   fList4 = fAnalysis->GetAcceptedDCAList();
170
171   fList5 = new TList();
172   fList5 = fAnalysis->GetEfficiencyQAList();
173 }
174
175 //________________________________________________________________________
176 void AliAnalysisTaskProtonsQA::Exec(Option_t *) {
177   // Main loop
178   // Called for each event
179   
180   if (!fESD) {
181     Printf("ERROR: fESD not available");
182     return;
183   }
184
185   if (!fMC) {
186     Printf("ERROR: Could not retrieve MC event");
187     return;
188   }
189   
190   AliStack* stack = fMC->Stack();
191   if (!stack) {
192     Printf("ERROR: Could not retrieve the stack");
193     return;
194   }
195   
196   fAnalysis->RunQAAnalysis(stack, fESD);
197   fAnalysis->RunMCAnalysis(stack);
198   fAnalysis->RunEfficiencyAnalysis(stack, fESD);
199
200   // Post output data.
201   PostData(0, fList0);
202   PostData(1, fList1);
203   PostData(2, fList2);
204   PostData(3, fList3);
205   PostData(4, fList4);
206   PostData(5, fList5);
207 }      
208
209 //________________________________________________________________________
210 void AliAnalysisTaskProtonsQA::Terminate(Option_t *) {
211   // Draw result to the screen
212   // Called once at the end of the query
213   
214   fList0 = dynamic_cast<TList*> (GetOutputData(0));
215   if (!fList0) {
216     Printf("ERROR: fList not available");
217     return;
218   }
219 }
220
221