]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/Correlations/DPhi/DiHadronPID/AliAODEventCutsDiHadronPID.cxx
DiHadronPID task update (Misha.Veldhoen@cern.ch)
[u/mrichter/AliRoot.git] / PWGCF / Correlations / DPhi / DiHadronPID / AliAODEventCutsDiHadronPID.cxx
1 // -------------------------------------------------------------------------
2 // Object managing event cuts, and holding QA histograms.
3 // -------------------------------------------------------------------------
4 // Author: Misha Veldhoen (misha.veldhoen@cern.ch)
5
6 #include <iostream>
7 #include "TList.h"
8 #include "TMath.h"
9 #include "TH1F.h"
10 #include "TH2F.h"
11 #include "TNamed.h"
12 #include "TIterator.h"
13 #include "AliAODEvent.h"
14 #include "AliAODHeader.h"
15 #include "AliAODVertex.h"
16 #include "AliCentrality.h"
17 #include "AliAnalysisManager.h"
18 #include "AliInputEventHandler.h"
19 #include "AliAODEventCutsDiHadronPID.h"
20
21 using namespace std;
22
23 ClassImp(AliAODEventCutsDiHadronPID);
24
25 // -------------------------------------------------------------------------
26 AliAODEventCutsDiHadronPID::AliAODEventCutsDiHadronPID():
27         TNamed(),
28         fIsPbPb(kTRUE),
29         fIsMC(kFALSE),
30         fTrigger(AliVEvent::kMB),
31         fMinCentrality(5.),
32         fMaxCentrality(0.),
33         fCentralityEstimator("V0M"),
34         fMaxVertexZ(10.),
35         fMinRefMult(0),
36         fTestTrigger(kFALSE),
37         fTestCentrality(kFALSE),
38         fTestVertexZ(kFALSE),
39         fTestMinRefMult(kFALSE),
40         fSelectedEventQAHistos(0x0),
41         fAllEventQAHistos(0x0),
42         fHistTrigger(0x0),
43         fHistRefMultiplicity(0x0),
44         fHistCentrality(0x0),
45         fHistCentralityQuality(0x0),
46         fHistVertexZ(0x0),
47         fDebug(0)
48
49 {
50
51         // 
52         // Default Constructor
53         //
54
55         cout<<"AliAODEventCutsDiHadronPID Default Constructor Called."<<endl;
56         if (fDebug) {cout << Form("File: %s, Line: %i, Function: %s",__FILE__,__LINE__,__func__) << endl;}
57
58 }
59
60 // -------------------------------------------------------------------------
61 AliAODEventCutsDiHadronPID::AliAODEventCutsDiHadronPID(const char* name):
62         TNamed(name,"AOD Event Cuts"),
63         fIsPbPb(kTRUE),
64         fIsMC(kFALSE),  
65         fTrigger(AliVEvent::kMB),
66         fMinCentrality(5.),
67         fMaxCentrality(0.),
68         fCentralityEstimator("V0M"),
69         fMaxVertexZ(10.),
70         fMinRefMult(0), 
71         fTestTrigger(kFALSE),
72         fTestCentrality(kFALSE),
73         fTestVertexZ(kFALSE),
74         fTestMinRefMult(kFALSE),        
75         fSelectedEventQAHistos(0x0),
76         fAllEventQAHistos(0x0),
77         fHistTrigger(0x0),
78         fHistRefMultiplicity(0x0),
79         fHistCentrality(0x0),
80         fHistCentralityQuality(0x0),
81         fHistVertexZ(0x0),
82         fDebug(0)
83         
84 {
85
86         //
87         // Named Constructor
88         //
89
90         cout<<"AliAODEventCutsDiHadronPID Named Constructor Called."<<endl;
91         if (fDebug) {cout << Form("File: %s, Line: %i, Function: %s",__FILE__,__LINE__,__func__) << endl;}
92
93 }
94
95 // -------------------------------------------------------------------------
96 AliAODEventCutsDiHadronPID::AliAODEventCutsDiHadronPID(const AliAODEventCutsDiHadronPID &source):
97         TNamed(source),
98         fTrigger(source.fTrigger),
99         fMinCentrality(source.fMinCentrality),
100         fMaxCentrality(source.fMaxCentrality),
101         fCentralityEstimator(source.fCentralityEstimator),
102         fMaxVertexZ(source.fMaxVertexZ),
103         fTestTrigger(source.fTestTrigger),
104         fTestCentrality(source.fTestCentrality),
105         fTestVertexZ(source.fTestVertexZ),
106         fSelectedEventQAHistos(0x0),
107         fAllEventQAHistos(0x0),
108         fHistTrigger(0x0),
109         fHistRefMultiplicity(0x0),
110         fHistCentrality(0x0),
111         fHistVertexZ(0x0)
112
113 {
114
115         //
116         // Copy Constructor
117         //
118
119         cout<<"AliAODEventCutsDiHadronPID Copy Constructor Called."<<endl;
120         if (fDebug) {cout << Form("File: %s, Line: %i, Function: %s",__FILE__,__LINE__,__func__) << endl;}
121
122         //if (source.fSelectedEventQAHistos) fSelectedEventQAHistos = (TList*)source.fHistTrigger->Clone();
123         //if (source.fHistTrigger) fHistTrigger = (TH1F*)source.fHistTrigger->Clone();
124         //if (source.fHistRefMultiplicity) fHistRefMultiplicity = (TH1F*)source.fHistRefMultiplicity->Clone();
125         //if (source.fHistCentrality) fHistCentrality = (TH1F*)source.fHistCentrality->Clone();
126         //if (source.fHistVertexZ) fHistVertexZ = (TH1F*)source.fHistVertexZ->Clone();
127
128 }
129
130 // -------------------------------------------------------------------------
131 AliAODEventCutsDiHadronPID::~AliAODEventCutsDiHadronPID() {
132
133         //
134         // Destructor
135         //
136
137         cout<<"AliAODEventCutsDiHadronPID Destructor Called."<<endl;
138         if (fDebug) {cout << Form("File: %s, Line: %i, Function: %s",__FILE__,__LINE__,__func__) << endl;}
139
140         if (fSelectedEventQAHistos) delete fSelectedEventQAHistos;
141         fSelectedEventQAHistos = 0x0;
142         if (fAllEventQAHistos) delete fAllEventQAHistos;
143         fAllEventQAHistos = 0x0;
144
145 }
146
147 // -------------------------------------------------------------------------
148 Long64_t AliAODEventCutsDiHadronPID::Merge(TCollection* list) {
149
150         //
151         // Merger. 
152         // 
153
154         cout<<"AliAODEventCutsDiHadronPID Merger Called."<<endl;
155         if (fDebug) {cout << Form("File: %s, Line: %i, Function: %s",__FILE__,__LINE__,__func__) << endl;}
156
157         if (!list) return 0;
158         if (list->IsEmpty()) return 1;
159
160         if (!fSelectedEventQAHistos||!fAllEventQAHistos) {
161                 cout<<"AliAODEventCutsDiHadronPID::Merge() - Warning, current object's histograms are missing... Generating."<<endl;
162                 CreateHistos();
163         }
164
165         TIterator* iter = list->MakeIterator();
166         TObject* obj;
167
168         // List of collections
169         TList collection_fSelectedEventQAHistos;
170         TList collection_fAllEventQAHistos;
171
172         Int_t count = 0;
173
174         while ((obj = iter->Next())) {
175         AliAODEventCutsDiHadronPID* entry = dynamic_cast<AliAODEventCutsDiHadronPID*> (obj);
176         if (entry == 0) continue;
177
178         // Check if the object to be merged really has the same name! (FIXME!)
179
180         // Getting the lists from obj.
181         TList* list_fSelectedEventQAHistos = entry->GetListOfSelectedEventQAHistos();
182         TList* list_fAllEventQAHistos = entry->GetListOfAllEventQAHistos();
183
184         // Adding the retrieved lists to the collection.
185         if (list_fSelectedEventQAHistos) collection_fSelectedEventQAHistos.Add(list_fSelectedEventQAHistos);
186         if (list_fAllEventQAHistos) collection_fAllEventQAHistos.Add(list_fAllEventQAHistos);
187
188         count++;
189     }
190
191     // Merging. Note that we require the original list to exist.
192     //  * Assume that if the collection happens to be empty, then nothing will happen.
193     //  * All other variables are taken from the original object.
194     if (fSelectedEventQAHistos) fSelectedEventQAHistos->Merge(&collection_fSelectedEventQAHistos);
195     if (fAllEventQAHistos) fAllEventQAHistos->Merge(&collection_fAllEventQAHistos);
196
197     delete iter;
198
199         return count+1;
200
201 }
202
203 // -------------------------------------------------------------------------
204 void AliAODEventCutsDiHadronPID::CreateHistos() {
205
206         cout<<"AliAODEventCutsDiHadronPID - Creating histograms"<<endl;
207         if (fDebug) {cout << Form("File: %s, Line: %i, Function: %s",__FILE__,__LINE__,__func__) << endl;}
208
209         // Create list of Event related QA histograms (selected events).
210         fSelectedEventQAHistos = new TList();
211         fSelectedEventQAHistos->SetName("SelectedEventQAHistos");
212         fSelectedEventQAHistos->SetOwner(kTRUE);
213
214         // The same, but for all events.
215         fAllEventQAHistos = new TList();
216         fAllEventQAHistos->SetName("AllEventQAHistos");
217         fAllEventQAHistos->SetOwner(kTRUE);
218
219         // Creating arrays of pointers to the QA histos.
220         fHistTrigger = new TH1F*[2];
221         fHistRefMultiplicity = new TH1F*[2];
222         fHistCentrality = new TH1F*[2];
223         fHistCentralityQuality = new TH1F*[2];
224         fHistVertexZ = new TH1F*[2];
225
226         const char* HistType[2] = {"Selected","All"};
227
228         for (Int_t iHistType = 0; iHistType < 2; iHistType++) {
229
230                 // Trigger Histogram.
231                 fHistTrigger[iHistType] = new TH1F(Form("fHistTrigger%s",HistType[iHistType]),"Trigger;;Count",2,-0.5,1.5);
232                 fHistTrigger[iHistType]->GetXaxis()->SetBinLabel(1,"kMB");
233                 fHistTrigger[iHistType]->GetXaxis()->SetBinLabel(2,"Other");
234                 if (iHistType == 0) fSelectedEventQAHistos->Add(fHistTrigger[iHistType]);
235                 else fAllEventQAHistos->Add(fHistTrigger[iHistType]);
236
237                 // Ref Multiplicity Histogram.
238                 fHistRefMultiplicity[iHistType] = new TH1F(Form("fHistRefMultiplicity%s",HistType[iHistType]),"Reference Multiplicity;N_{tracks};Count",100,0.,10000.);
239                 if (iHistType == 0) fSelectedEventQAHistos->Add(fHistRefMultiplicity[iHistType]);
240                 else fAllEventQAHistos->Add(fHistRefMultiplicity[iHistType]);
241
242                 // Centrality Histogram.
243                 fHistCentrality[iHistType] = new TH1F(Form("fHistCentrality%s",HistType[iHistType]),"Centrality;Centrality;Count",20,0,100);
244                 if (iHistType == 0) fSelectedEventQAHistos->Add(fHistCentrality[iHistType]);
245                 else fAllEventQAHistos->Add(fHistCentrality[iHistType]);
246
247                 // Centrality Quality.
248                 fHistCentralityQuality[iHistType] = new TH1F(Form("fHistCentralityQuality%s",HistType[iHistType]),"Centrality Quality;Quality;Count",2,-0.5,1.5);
249                 fHistCentralityQuality[iHistType]->GetXaxis()->SetBinLabel(1,"0");
250                 fHistCentralityQuality[iHistType]->GetXaxis()->SetBinLabel(2,"Other");
251                 if (iHistType == 0) fSelectedEventQAHistos->Add(fHistCentralityQuality[iHistType]);
252                 else fAllEventQAHistos->Add(fHistCentralityQuality[iHistType]);
253
254                 // VertexZ Histogram.
255                 fHistVertexZ[iHistType] = new TH1F(Form("fHistVertexZ%s",HistType[iHistType]),"VertexZ;z (cm);Count",60,-15.,15.);
256                 if (iHistType == 0) fSelectedEventQAHistos->Add(fHistVertexZ[iHistType]);
257                 else fAllEventQAHistos->Add(fHistVertexZ[iHistType]);
258
259         }
260
261 }
262
263 // -------------------------------------------------------------------------
264 Bool_t AliAODEventCutsDiHadronPID::IsSelected(AliAODEvent* event) {
265         
266         if (fDebug) {cout << Form("File: %s, Line: %i, Function: %s",__FILE__,__LINE__,__func__) << endl;}
267
268         if (!event) return kFALSE;
269
270         if (!fAllEventQAHistos||!fSelectedEventQAHistos) {cout<<"AliAODEventCutsDiHadronPID - Histograms were not created, you should have called CreateHistos()..."<<endl;}
271
272         // Input the event handler.
273         AliInputEventHandler* InputHandler = (AliInputEventHandler*)((AliAnalysisManager::GetAnalysisManager())->GetInputEventHandler());
274         if (!InputHandler) return kFALSE;
275
276         Bool_t select = kTRUE;
277
278         // Test Trigger.
279         UInt_t trigger = InputHandler->IsEventSelected();
280         Int_t triggerselect = 0; // 0 = selected.
281     if (fTestTrigger) {
282                 if (!(trigger & fTrigger)) {
283                         select = kFALSE;
284                         triggerselect = 1;      // 1 = not selected.
285                 }
286         }       
287
288         AliCentrality* CurrentCentrality = 0x0;
289         Int_t CurrentCentralityQuality = -999;
290         Float_t percentile = -999.;
291
292         if (fIsPbPb) {
293
294                 // Get the centrality object.
295             CurrentCentrality = event->GetCentrality();
296             if (!CurrentCentrality) select = kFALSE;
297
298             // Check the quality of the centrality estimation.
299             // If 0 then quality is OK, c.f. TOF/PbPb276/macros/TOFmatchEff.C
300             CurrentCentralityQuality = CurrentCentrality->GetQuality();
301             //cout<<"Centrality: "<<CurrentCentrality->GetCentralityPercentile(fCentralityEstimator.Data())<<" Quality: "<<CurrentCentrality->GetQuality()<<endl;
302             if (CurrentCentralityQuality) select = kFALSE;
303
304                 // Test Centrality.
305             percentile = CurrentCentrality->GetCentralityPercentile(fCentralityEstimator.Data());
306                 if (fTestCentrality) {
307                 if ((percentile < fMaxCentrality)||(percentile > fMinCentrality)) select = kFALSE;
308                 }
309
310         }
311
312         // Get the primary vertex.
313         AliAODVertex* CurrentPrimaryVertex = event->GetPrimaryVertex();
314     if (!CurrentPrimaryVertex) select = kFALSE;
315
316         // Test Vertex Z.
317     Double_t vtxz = CurrentPrimaryVertex->GetZ();
318         if (fTestVertexZ) {
319         if (TMath::Abs(vtxz) > fMaxVertexZ) select = kFALSE; 
320         }
321
322         // Get the event header.
323         AliAODHeader* CurrentHeader = event->GetHeader();
324
325         // Test minimum reference multiplicity.
326         Int_t CurrentRefMultiplicity = CurrentHeader->GetRefMultiplicity();
327         if (fTestMinRefMult) {
328                 if (CurrentRefMultiplicity < fMinRefMult) select = kFALSE;
329         }
330
331         // Fill the histograms for selected events.
332         if (select) {
333                 fHistTrigger[0]->Fill(triggerselect);
334                 fHistRefMultiplicity[0]->Fill(CurrentHeader->GetRefMultiplicity());
335                 if (fIsPbPb) fHistCentrality[0]->Fill(percentile);
336                 if (fIsPbPb) fHistCentralityQuality[0]->Fill(CurrentCentralityQuality);
337                 fHistVertexZ[0]->Fill(vtxz);
338         }
339
340         // Fill the histograms for all events.
341         fHistTrigger[1]->Fill(triggerselect);
342         fHistRefMultiplicity[1]->Fill(CurrentHeader->GetRefMultiplicity());
343         if (fIsPbPb) fHistCentrality[1]->Fill(percentile);
344         if (fIsPbPb) fHistCentralityQuality[1]->Fill(CurrentCentralityQuality);
345         fHistVertexZ[1]->Fill(vtxz);
346
347         cout<<"Event Selected: "<<select<<endl;
348
349         return select;
350
351 }
352
353 // -------------------------------------------------------------------------
354 void AliAODEventCutsDiHadronPID::PrintCuts() {
355
356         // NOT IMPLEMENTED.
357         if (fDebug) {cout << Form("File: %s, Line: %i, Function: %s",__FILE__,__LINE__,__func__) << endl;}
358
359         return;
360
361 }