]>
Commit | Line | Data |
---|---|---|
563113d0 | 1 | /************************************************************************** |
2 | * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. * | |
3 | * * | |
4 | * Author: The ALICE Off-line Project. * | |
5 | * Contributors are mentioned in the code where appropriate. * | |
6 | * * | |
7 | * Permission to use, copy, modify and distribute this software and its * | |
8 | * documentation strictly for non-commercial purposes is hereby granted * | |
9 | * without fee, provided that the above copyright notice appears in all * | |
10 | * copies and that both the copyright notice and this permission notice * | |
11 | * appear in the supporting documentation. The authors make no claims * | |
12 | * about the suitability of this software for any purpose. It is * | |
13 | * provided "as is" without express or implied warranty. * | |
14 | **************************************************************************/ | |
15 | ||
16 | // The class AliCFTrackKineCuts is designed to select both generated | |
17 | // and reconstructed tracks of a given range in momentum space, | |
18 | // electric charge and azimuthal emission angle phi | |
19 | // and to provide corresponding QA histograms. | |
20 | // This class inherits from the Analysis' Framework abstract base class | |
21 | // AliAnalysisCuts and is a part of the Correction Framework. | |
22 | // This class acts on single, generated and reconstructed tracks, it is | |
23 | // applicable on ESD and AOD data. | |
24 | // It mainly consists of a IsSelected function that returns a boolean. | |
25 | // This function checks whether the considered track passes a set of cuts: | |
26 | // - total momentum | |
27 | // - pt | |
28 | // - px | |
29 | // - py | |
30 | // - pz | |
31 | // - eta | |
32 | // - rapidity | |
33 | // - phi | |
34 | // - charge | |
35 | // - is charged | |
36 | // | |
37 | // The cut values for these cuts are set with the corresponding set functions. | |
38 | // All cut classes provided by the correction framework are supposed to be | |
39 | // added in the Analysis Framwork's class AliAnalysisFilter and applied by | |
40 | // the filter via a loop. | |
41 | // | |
42 | // author: I. Kraus (Ingrid.Kraus@cern.ch) | |
43 | // idea taken form | |
44 | // AliESDtrackCuts writte by Jan Fiete Grosse-Oetringhaus and | |
45 | // AliRsnDaughterCut class written by A. Pulvirenti. | |
46 | ||
47 | #include <TCanvas.h> | |
48 | #include <TDirectory.h> | |
49 | #include <TBits.h> | |
50 | #include <TH2.h> | |
51 | ||
52 | #include <AliVParticle.h> | |
53 | #include <AliLog.h> | |
54 | #include "AliCFTrackKineCuts.h" | |
55 | ||
56 | ClassImp(AliCFTrackKineCuts) | |
57 | ||
58 | //__________________________________________________________________________________ | |
59 | AliCFTrackKineCuts::AliCFTrackKineCuts() : | |
60 | AliCFCutBase(), | |
61 | fMomentumMin(0), | |
62 | fMomentumMax(0), | |
63 | fPtMin(0), | |
64 | fPtMax(0), | |
65 | fPxMin(0), | |
66 | fPxMax(0), | |
67 | fPyMin(0), | |
68 | fPyMax(0), | |
69 | fPzMin(0), | |
70 | fPzMax(0), | |
71 | fEtaMin(0), | |
72 | fEtaMax(0), | |
73 | fRapidityMin(0), | |
74 | fRapidityMax(0), | |
75 | fPhiMin(0), | |
76 | fPhiMax(0), | |
77 | fCharge(0), | |
78 | fRequireIsCharged(0), | |
79 | fhCutStatistics(0), | |
80 | fhCutCorrelation(0), | |
81 | fBitmap(0x0), | |
82 | fhNBinsMomentum(0), | |
83 | fhNBinsPt(0), | |
84 | fhNBinsPx(0), | |
85 | fhNBinsPy(0), | |
86 | fhNBinsPz(0), | |
87 | fhNBinsEta(0), | |
88 | fhNBinsRapidity(0), | |
89 | fhNBinsPhi(0), | |
90 | fhNBinsCharge(0), | |
91 | fhBinLimMomentum(0x0), | |
92 | fhBinLimPt(0x0), | |
93 | fhBinLimPx(0x0), | |
94 | fhBinLimPy(0x0), | |
95 | fhBinLimPz(0x0), | |
96 | fhBinLimEta(0x0), | |
97 | fhBinLimRapidity(0x0), | |
98 | fhBinLimPhi(0x0), | |
99 | fhBinLimCharge(0x0) | |
100 | ||
101 | { | |
102 | // | |
103 | // Default constructor | |
104 | // | |
105 | fBitmap=new TBits(0); | |
106 | Initialise(); | |
107 | } | |
108 | //__________________________________________________________________________________ | |
109 | AliCFTrackKineCuts::AliCFTrackKineCuts(Char_t* name, Char_t* title) : | |
110 | AliCFCutBase(name,title), | |
111 | fMomentumMin(0), | |
112 | fMomentumMax(0), | |
113 | fPtMin(0), | |
114 | fPtMax(0), | |
115 | fPxMin(0), | |
116 | fPxMax(0), | |
117 | fPyMin(0), | |
118 | fPyMax(0), | |
119 | fPzMin(0), | |
120 | fPzMax(0), | |
121 | fEtaMin(0), | |
122 | fEtaMax(0), | |
123 | fRapidityMin(0), | |
124 | fRapidityMax(0), | |
125 | fPhiMin(0), | |
126 | fPhiMax(0), | |
127 | fCharge(0), | |
128 | fRequireIsCharged(0), | |
129 | fhCutStatistics(0), | |
130 | fhCutCorrelation(0), | |
131 | fBitmap(0x0), | |
132 | fhNBinsMomentum(0), | |
133 | fhNBinsPt(0), | |
134 | fhNBinsPx(0), | |
135 | fhNBinsPy(0), | |
136 | fhNBinsPz(0), | |
137 | fhNBinsEta(0), | |
138 | fhNBinsRapidity(0), | |
139 | fhNBinsPhi(0), | |
140 | fhNBinsCharge(0), | |
141 | fhBinLimMomentum(0x0), | |
142 | fhBinLimPt(0x0), | |
143 | fhBinLimPx(0x0), | |
144 | fhBinLimPy(0x0), | |
145 | fhBinLimPz(0x0), | |
146 | fhBinLimEta(0x0), | |
147 | fhBinLimRapidity(0x0), | |
148 | fhBinLimPhi(0x0), | |
149 | fhBinLimCharge(0x0) | |
150 | ||
151 | { | |
152 | // | |
153 | // Constructor | |
154 | // | |
155 | fBitmap=new TBits(0); | |
156 | Initialise(); | |
157 | } | |
158 | //__________________________________________________________________________________ | |
159 | AliCFTrackKineCuts::AliCFTrackKineCuts(const AliCFTrackKineCuts& c) : | |
160 | AliCFCutBase(c), | |
161 | fMomentumMin(c.fMomentumMin), | |
162 | fMomentumMax(c.fMomentumMax), | |
163 | fPtMin(c.fPtMin), | |
164 | fPtMax(c.fPtMax), | |
165 | fPxMin(c.fPxMin), | |
166 | fPxMax(c.fPxMax), | |
167 | fPyMin(c.fPyMin), | |
168 | fPyMax(c.fPyMax), | |
169 | fPzMin(c.fPzMin), | |
170 | fPzMax(c.fPzMax), | |
171 | fEtaMin(c.fEtaMin), | |
172 | fEtaMax(c.fEtaMax), | |
173 | fRapidityMin(c.fRapidityMin), | |
174 | fRapidityMax(c.fRapidityMax), | |
175 | fPhiMin(c.fPhiMin), | |
176 | fPhiMax(c.fPhiMax), | |
177 | fCharge(c.fCharge), | |
178 | fRequireIsCharged(c.fRequireIsCharged), | |
179 | fhCutStatistics(c.fhCutStatistics), | |
180 | fhCutCorrelation(c.fhCutCorrelation), | |
181 | fBitmap(c.fBitmap), | |
182 | fhNBinsMomentum(c.fhNBinsMomentum), | |
183 | fhNBinsPt(c.fhNBinsPt), | |
184 | fhNBinsPx(c.fhNBinsPx), | |
185 | fhNBinsPy(c.fhNBinsPy), | |
186 | fhNBinsPz(c.fhNBinsPz), | |
187 | fhNBinsEta(c.fhNBinsEta), | |
188 | fhNBinsRapidity(c.fhNBinsRapidity), | |
189 | fhNBinsPhi(c.fhNBinsPhi), | |
190 | fhNBinsCharge(c.fhNBinsCharge), | |
191 | fhBinLimMomentum(c.fhBinLimMomentum), | |
192 | fhBinLimPt(c.fhBinLimPt), | |
193 | fhBinLimPx(c.fhBinLimPx), | |
194 | fhBinLimPy(c.fhBinLimPy), | |
195 | fhBinLimPz(c.fhBinLimPz), | |
196 | fhBinLimEta(c.fhBinLimEta), | |
197 | fhBinLimRapidity(c.fhBinLimRapidity), | |
198 | fhBinLimPhi(c.fhBinLimPhi), | |
199 | fhBinLimCharge(c.fhBinLimCharge) | |
200 | ||
201 | { | |
202 | // | |
203 | // copy constructor | |
204 | // | |
205 | ((AliCFTrackKineCuts &) c).Copy(*this); | |
206 | } | |
207 | //__________________________________________________________________________________ | |
208 | AliCFTrackKineCuts& AliCFTrackKineCuts::operator=(const AliCFTrackKineCuts& c) | |
209 | { | |
210 | // | |
211 | // Assignment operator | |
212 | // | |
213 | if (this != &c) { | |
214 | AliCFCutBase::operator=(c) ; | |
215 | fMomentumMin = c.fMomentumMin ; | |
216 | fMomentumMax = c.fMomentumMax ; | |
217 | fPtMin = c.fPtMin ; | |
218 | fPtMax = c.fPtMax ; | |
219 | fPxMin = c.fPxMin ; | |
220 | fPxMax = c.fPxMax ; | |
221 | fPyMin = c.fPyMin ; | |
222 | fPyMax = c.fPyMax ; | |
223 | fPzMin = c.fPzMin ; | |
224 | fPzMax = c.fPzMax ; | |
225 | fEtaMin = c.fEtaMin ; | |
226 | fEtaMax = c.fEtaMax ; | |
227 | fRapidityMin = c.fRapidityMin ; | |
228 | fRapidityMax = c.fRapidityMax ; | |
229 | fPhiMin = c.fPhiMin ; | |
230 | fPhiMax = c.fPhiMax ; | |
231 | fCharge = c.fCharge ; | |
232 | fRequireIsCharged = c.fRequireIsCharged ; | |
233 | fhCutStatistics = c.fhCutStatistics ; | |
234 | fhCutCorrelation = c.fhCutCorrelation ; | |
235 | fBitmap = c.fBitmap; | |
236 | fhNBinsMomentum = c.fhNBinsMomentum; | |
237 | fhNBinsPt = c.fhNBinsPt; | |
238 | fhNBinsPx = c.fhNBinsPx; | |
239 | fhNBinsPy = c.fhNBinsPy; | |
240 | fhNBinsPz = c.fhNBinsPz; | |
241 | fhNBinsEta = c.fhNBinsEta; | |
242 | fhNBinsRapidity = c.fhNBinsRapidity; | |
243 | fhNBinsPhi = c.fhNBinsPhi; | |
244 | fhNBinsCharge = c.fhNBinsCharge; | |
245 | fhBinLimMomentum = c.fhBinLimMomentum; | |
246 | fhBinLimPt = c.fhBinLimPt; | |
247 | fhBinLimPx = c.fhBinLimPx; | |
248 | fhBinLimPy = c.fhBinLimPy; | |
249 | fhBinLimPz = c.fhBinLimPz; | |
250 | fhBinLimEta = c.fhBinLimEta; | |
251 | fhBinLimRapidity = c.fhBinLimRapidity; | |
252 | fhBinLimPhi = c.fhBinLimPhi; | |
253 | fhBinLimCharge = c.fhBinLimCharge; | |
254 | ||
255 | for (Int_t i=0; i<c.kNHist; i++){ | |
256 | for (Int_t j=0; j<c.kNStepQA; j++){ | |
257 | if(c.fhQA[i][j]) fhQA[i][j] = (TH1F*)c.fhQA[i][j]->Clone(); | |
258 | } | |
259 | } | |
260 | ||
261 | ((AliCFTrackKineCuts &) c).Copy(*this); | |
262 | } | |
263 | return *this ; | |
264 | } | |
265 | //__________________________________________________________________________________ | |
266 | AliCFTrackKineCuts::~AliCFTrackKineCuts() | |
267 | { | |
268 | // | |
269 | // destructor | |
270 | // | |
271 | if (fhCutStatistics) delete fhCutStatistics; | |
272 | if (fhCutCorrelation) delete fhCutCorrelation; | |
273 | ||
274 | for (Int_t i=0; i<kNHist; i++){ | |
275 | for (Int_t j=0; j<kNStepQA; j++){ | |
276 | if(fhQA[i][j]) delete fhQA[i][j]; | |
277 | } | |
278 | } | |
279 | ||
280 | if(fBitmap) delete fBitmap; | |
281 | ||
282 | if(fhBinLimMomentum) delete fhBinLimMomentum; | |
283 | if(fhBinLimPt) delete fhBinLimPt; | |
284 | if(fhBinLimPx) delete fhBinLimPx; | |
285 | if(fhBinLimPy) delete fhBinLimPy; | |
286 | if(fhBinLimPz) delete fhBinLimPz; | |
287 | if(fhBinLimEta) delete fhBinLimEta; | |
288 | if(fhBinLimRapidity) delete fhBinLimRapidity; | |
289 | if(fhBinLimPhi) delete fhBinLimPhi; | |
290 | if(fhBinLimCharge) delete fhBinLimCharge; | |
291 | } | |
292 | //__________________________________________________________________________________ | |
293 | void AliCFTrackKineCuts::Initialise() | |
294 | { | |
295 | // | |
296 | // sets everything to zero | |
297 | // | |
298 | fMomentumMin = 0; | |
299 | fMomentumMax = 0; | |
300 | fPtMin = 0; | |
301 | fPtMax = 0; | |
302 | fPxMin = 0; | |
303 | fPxMax = 0; | |
304 | fPyMin = 0; | |
305 | fPyMax = 0; | |
306 | fPzMin = 0; | |
307 | fPzMax = 0; | |
308 | fEtaMin = 0; | |
309 | fEtaMax = 0; | |
310 | fRapidityMin = 0; | |
311 | fRapidityMax = 0; | |
312 | fPhiMin = 0; | |
313 | fPhiMax = 0; | |
314 | fCharge = 0; | |
315 | fRequireIsCharged = 0; | |
316 | ||
317 | SetMomentumRange(); | |
318 | SetPtRange(); | |
319 | SetPxRange(); | |
320 | SetPyRange(); | |
321 | SetPzRange(); | |
322 | SetEtaRange(); | |
323 | SetRapidityRange(); | |
324 | SetPhiRange(); | |
325 | SetChargeRec(); | |
326 | SetChargeMC(); | |
327 | SetRequireIsCharged(); | |
328 | ||
329 | for (Int_t i=0; i<kNHist; i++){ | |
330 | for (Int_t j=0; j<kNStepQA; j++) { | |
331 | fhQA[i][j] = 0x0; | |
332 | } | |
333 | } | |
334 | ||
335 | fhCutStatistics = 0; | |
336 | fhCutCorrelation = 0; | |
337 | ||
338 | //set default bining for QA histograms | |
339 | SetHistogramBins(kCutP,200,0.,20.); | |
340 | SetHistogramBins(kCutPt,200,0.,20.); | |
341 | SetHistogramBins(kCutPx,400,-20.,20.); | |
342 | SetHistogramBins(kCutPy,400,-20.,20.); | |
343 | SetHistogramBins(kCutPz,400,-20.,20.); | |
344 | SetHistogramBins(kCutRapidity,200,-10.,10.); | |
345 | SetHistogramBins(kCutEta,200,-10.,10.); | |
346 | SetHistogramBins(kCutPhi,38,-0.6,7.); | |
347 | SetHistogramBins(kCutCharge,21,-10.5,10.5); | |
348 | ||
349 | } | |
350 | //__________________________________________________________________________________ | |
351 | void AliCFTrackKineCuts::Copy(TObject &c) const | |
352 | { | |
353 | // | |
354 | // Copy function | |
355 | // | |
356 | AliCFTrackKineCuts& target = (AliCFTrackKineCuts &) c; | |
357 | ||
358 | target.Initialise(); | |
359 | ||
360 | if (fhCutStatistics) target.fhCutStatistics = (TH1F*) fhCutStatistics->Clone(); | |
361 | if (fhCutCorrelation) target.fhCutCorrelation = (TH2F*) fhCutCorrelation->Clone(); | |
362 | ||
363 | for (Int_t i=0; i<kNHist; i++){ | |
364 | for (Int_t j=0; j<kNStepQA; j++){ | |
365 | if(fhQA[i][j]) target.fhQA[i][j] = (TH1F*)fhQA[i][j]->Clone(); | |
366 | } | |
367 | } | |
368 | ||
369 | TNamed::Copy(c); | |
370 | } | |
371 | //__________________________________________________________________________________ | |
372 | void AliCFTrackKineCuts::GetBitMap(TObject* obj, TBits *bitmap) { | |
373 | // | |
374 | // retrieve the pointer to the bitmap | |
375 | // | |
376 | bitmap = SelectionBitMap(obj); | |
377 | } | |
378 | //__________________________________________________________________________________ | |
379 | TBits* AliCFTrackKineCuts::SelectionBitMap(TObject* obj) { | |
380 | // | |
381 | // test if the track passes the single cuts | |
382 | // and store the information in a bitmap | |
383 | // | |
384 | ||
385 | // bitmap stores the decision of each single cut | |
386 | for(Int_t i=0; i<kNCuts; i++)fBitmap->SetBitNumber(i,kFALSE); | |
387 | // cast TObject into VParticle | |
388 | AliVParticle* particle = dynamic_cast<AliVParticle *>(obj); | |
389 | if ( !particle ) return fBitmap ; | |
390 | ||
391 | ||
392 | for(Int_t i=0; i<kNCuts; i++)fBitmap->SetBitNumber(i,kTRUE); | |
393 | ||
394 | Int_t iCutBit = 0; | |
395 | if((particle->P() < fMomentumMin) || (particle->P() > fMomentumMax)) | |
396 | fBitmap->SetBitNumber(iCutBit,kFALSE); | |
397 | iCutBit++; | |
398 | if ((particle->Pt() < fPtMin) || (particle->Pt() > fPtMax)) | |
399 | fBitmap->SetBitNumber(iCutBit,kFALSE); | |
400 | iCutBit++; | |
401 | if ((particle->Px() < fPxMin) || (particle->Px() > fPxMax)) | |
402 | fBitmap->SetBitNumber(iCutBit,kFALSE); | |
403 | iCutBit++; | |
404 | if ((particle->Py() < fPyMin) || (particle->Py() > fPyMax)) | |
405 | fBitmap->SetBitNumber(iCutBit,kFALSE); | |
406 | iCutBit++; | |
407 | if ((particle->Pz() < fPzMin) || (particle->Pz() > fPzMax)) | |
408 | fBitmap->SetBitNumber(iCutBit,kFALSE); | |
409 | iCutBit++; | |
410 | if ((particle->Eta() < fEtaMin) || (particle->Eta() > fEtaMax)) | |
411 | fBitmap->SetBitNumber(iCutBit,kFALSE); | |
412 | iCutBit++; | |
413 | if ((particle->Y() < fRapidityMin) || (particle->Y() > fRapidityMax)) | |
414 | fBitmap->SetBitNumber(iCutBit,kFALSE); | |
415 | iCutBit++; | |
416 | if ((particle->Phi() < fPhiMin) || (particle->Phi() > fPhiMax)) | |
417 | fBitmap->SetBitNumber(iCutBit,kFALSE); | |
418 | iCutBit++; | |
419 | if (fCharge < 10 && particle->Charge() != fCharge) | |
420 | fBitmap->SetBitNumber(iCutBit,kFALSE); | |
421 | iCutBit++; | |
422 | if (fRequireIsCharged && particle->Charge()==0) | |
423 | fBitmap->SetBitNumber(iCutBit,kFALSE); | |
424 | ||
425 | return fBitmap; | |
426 | } | |
427 | //__________________________________________________________________________________ | |
428 | Bool_t AliCFTrackKineCuts::IsSelected(TObject* obj) { | |
429 | // | |
430 | // loops over decisions of single cuts and returns if the track is accepted | |
431 | // | |
432 | TBits* bitmap = SelectionBitMap(obj); | |
433 | ||
434 | Bool_t isSelected = kTRUE; | |
435 | ||
436 | for (UInt_t icut=0; icut<bitmap->GetNbits();icut++) | |
437 | if(!bitmap->TestBitNumber(icut)) isSelected = kFALSE; | |
438 | ||
439 | return isSelected; | |
440 | } | |
441 | //__________________________________________________________________________________ | |
442 | void AliCFTrackKineCuts::Init() { | |
443 | // | |
444 | // initialises all histograms and the TList which holds the histograms | |
445 | // | |
446 | if(fIsQAOn) | |
447 | DefineHistograms(); | |
448 | } | |
449 | //__________________________________________________________________________________ | |
450 | void AliCFTrackKineCuts::SetHistogramBins(Int_t index, Int_t nbins, Double_t *bins) | |
451 | { | |
452 | // | |
453 | // variable bin size | |
454 | // | |
455 | if(!fIsQAOn) return; | |
456 | ||
457 | switch(index){ | |
458 | case kCutP: | |
459 | fhNBinsMomentum=nbins; | |
460 | fhBinLimMomentum=new Double_t[nbins+1]; | |
461 | for(Int_t i=0;i<nbins+1;i++)fhBinLimMomentum[i]=bins[i]; | |
462 | break; | |
463 | ||
464 | case kCutPt: | |
465 | fhNBinsPt=nbins; | |
466 | fhBinLimPt=new Double_t[nbins+1]; | |
467 | for(Int_t i=0;i<nbins+1;i++)fhBinLimPt[i]=bins[i]; | |
468 | break; | |
469 | ||
470 | case kCutPx: | |
471 | fhNBinsPx=nbins; | |
472 | fhBinLimPx=new Double_t[nbins+1]; | |
473 | for(Int_t i=0;i<nbins+1;i++)fhBinLimPx[i]=bins[i]; | |
474 | break; | |
475 | ||
476 | case kCutPy: | |
477 | fhNBinsPy=nbins; | |
478 | fhBinLimPy=new Double_t[nbins+1]; | |
479 | for(Int_t i=0;i<nbins+1;i++)fhBinLimPy[i]=bins[i]; | |
480 | break; | |
481 | ||
482 | case kCutPz: | |
483 | fhNBinsPz=nbins; | |
484 | fhBinLimPz=new Double_t[nbins+1]; | |
485 | for(Int_t i=0;i<nbins+1;i++)fhBinLimPz[i]=bins[i]; | |
486 | break; | |
487 | ||
488 | case kCutRapidity: | |
489 | fhNBinsRapidity=nbins; | |
490 | fhBinLimRapidity=new Double_t[nbins+1]; | |
491 | for(Int_t i=0;i<nbins+1;i++)fhBinLimRapidity[i]=bins[i]; | |
492 | break; | |
493 | ||
494 | case kCutEta: | |
495 | fhNBinsEta=nbins; | |
496 | fhBinLimEta=new Double_t[nbins+1]; | |
497 | for(Int_t i=0;i<nbins+1;i++)fhBinLimEta[i]=bins[i]; | |
498 | break; | |
499 | ||
500 | case kCutPhi: | |
501 | fhNBinsPhi=nbins; | |
502 | fhBinLimPhi=new Double_t[nbins+1]; | |
503 | for(Int_t i=0;i<nbins+1;i++)fhBinLimPhi[i]=bins[i]; | |
504 | break; | |
505 | ||
506 | case kCutCharge: | |
507 | fhNBinsCharge=nbins; | |
508 | fhBinLimCharge=new Double_t[nbins+1]; | |
509 | for(Int_t i=0;i<nbins+1;i++)fhBinLimCharge[i]=bins[i]; | |
510 | break; | |
511 | } | |
512 | } | |
513 | //__________________________________________________________________________________ | |
514 | void AliCFTrackKineCuts::SetHistogramBins(Int_t index, Int_t nbins, Double_t xmin, Double_t xmax) | |
515 | { | |
516 | // | |
517 | // fixed bin size | |
518 | // | |
519 | if(!fIsQAOn) return; | |
520 | ||
521 | switch(index){ | |
522 | case kCutP: | |
523 | fhNBinsMomentum=nbins; | |
524 | fhBinLimMomentum=new Double_t[nbins+1]; | |
525 | for(Int_t i=0;i<nbins+1;i++)fhBinLimMomentum[i]=xmin+i*(xmax-xmin)/Double_t(nbins); | |
526 | break; | |
527 | ||
528 | case kCutPt: | |
529 | fhNBinsPt=nbins; | |
530 | fhBinLimPt=new Double_t[nbins+1]; | |
531 | for(Int_t i=0;i<nbins+1;i++)fhBinLimPt[i]=xmin+i*(xmax-xmin)/Double_t(nbins); | |
532 | break; | |
533 | ||
534 | case kCutPx: | |
535 | fhNBinsPx=nbins; | |
536 | fhBinLimPx=new Double_t[nbins+1]; | |
537 | for(Int_t i=0;i<nbins+1;i++)fhBinLimPx[i]=xmin+i*(xmax-xmin)/Double_t(nbins); | |
538 | break; | |
539 | ||
540 | case kCutPy: | |
541 | fhNBinsPy=nbins; | |
542 | fhBinLimPy=new Double_t[nbins+1]; | |
543 | for(Int_t i=0;i<nbins+1;i++)fhBinLimPy[i]=xmin+i*(xmax-xmin)/Double_t(nbins); | |
544 | break; | |
545 | ||
546 | case kCutPz: | |
547 | fhNBinsPz=nbins; | |
548 | fhBinLimPz=new Double_t[nbins+1]; | |
549 | for(Int_t i=0;i<nbins+1;i++)fhBinLimPz[i]=xmin+i*(xmax-xmin)/Double_t(nbins); | |
550 | break; | |
551 | ||
552 | case kCutRapidity: | |
553 | fhNBinsRapidity=nbins; | |
554 | fhBinLimRapidity=new Double_t[nbins+1]; | |
555 | for(Int_t i=0;i<nbins+1;i++)fhBinLimRapidity[i]=xmin+i*(xmax-xmin)/Double_t(nbins); | |
556 | break; | |
557 | ||
558 | case kCutEta: | |
559 | fhNBinsEta=nbins; | |
560 | fhBinLimEta=new Double_t[nbins+1]; | |
561 | for(Int_t i=0;i<nbins+1;i++)fhBinLimEta[i]=xmin+i*(xmax-xmin)/Double_t(nbins); | |
562 | break; | |
563 | ||
564 | case kCutPhi: | |
565 | fhNBinsPhi=nbins; | |
566 | fhBinLimPhi=new Double_t[nbins+1]; | |
567 | for(Int_t i=0;i<nbins+1;i++)fhBinLimPhi[i]=xmin+i*(xmax-xmin)/Double_t(nbins); | |
568 | break; | |
569 | ||
570 | case kCutCharge: | |
571 | fhNBinsCharge=nbins; | |
572 | fhBinLimCharge=new Double_t[nbins+1]; | |
573 | for(Int_t i=0;i<nbins+1;i++)fhBinLimCharge[i]=xmin+i*(xmax-xmin)/Double_t(nbins); | |
574 | break; | |
575 | } | |
576 | } | |
577 | //__________________________________________________________________________________ | |
578 | void AliCFTrackKineCuts::DefineHistograms() { | |
579 | // | |
580 | // histograms for cut variables, cut statistics and cut correlations | |
581 | // | |
582 | ||
583 | Int_t color = 2; | |
584 | ||
585 | // book cut statistics and cut correlation histograms | |
586 | fhCutStatistics = new TH1F(Form("%s_cut_statistics",GetName()), Form("%s cut statistics",GetName()), kNCuts,0.5,kNCuts+0.5); | |
587 | fhCutStatistics->SetLineWidth(2); | |
588 | fhCutStatistics->GetXaxis()->SetBinLabel(1,"p"); | |
589 | fhCutStatistics->GetXaxis()->SetBinLabel(2,"pt"); | |
590 | fhCutStatistics->GetXaxis()->SetBinLabel(3,"px"); | |
591 | fhCutStatistics->GetXaxis()->SetBinLabel(4,"py"); | |
592 | fhCutStatistics->GetXaxis()->SetBinLabel(5,"pz"); | |
593 | fhCutStatistics->GetXaxis()->SetBinLabel(6,"eta"); | |
594 | fhCutStatistics->GetXaxis()->SetBinLabel(7,"y"); | |
595 | fhCutStatistics->GetXaxis()->SetBinLabel(8,"phi"); | |
596 | fhCutStatistics->GetXaxis()->SetBinLabel(9,"charge"); | |
597 | fhCutStatistics->GetXaxis()->SetBinLabel(10,"is charged"); | |
598 | ||
599 | fhCutCorrelation = new TH2F(Form("%s_cut_correlation",GetName()), Form("%s cut correlation",GetName()), kNCuts,0.5,kNCuts+0.5,kNCuts,0.5,kNCuts+0.5); | |
600 | fhCutCorrelation->SetLineWidth(2); | |
601 | fhCutCorrelation->GetXaxis()->SetBinLabel(1,"p"); | |
602 | fhCutCorrelation->GetXaxis()->SetBinLabel(2,"pt"); | |
603 | fhCutCorrelation->GetXaxis()->SetBinLabel(3,"px"); | |
604 | fhCutCorrelation->GetXaxis()->SetBinLabel(4,"py"); | |
605 | fhCutCorrelation->GetXaxis()->SetBinLabel(5,"pz"); | |
606 | fhCutCorrelation->GetXaxis()->SetBinLabel(6,"eta"); | |
607 | fhCutCorrelation->GetXaxis()->SetBinLabel(7,"y"); | |
608 | fhCutCorrelation->GetXaxis()->SetBinLabel(8,"phi"); | |
609 | fhCutCorrelation->GetXaxis()->SetBinLabel(9,"charge"); | |
610 | fhCutCorrelation->GetXaxis()->SetBinLabel(10,"is charged"); | |
611 | ||
612 | fhCutCorrelation->GetYaxis()->SetBinLabel(1,"p"); | |
613 | fhCutCorrelation->GetYaxis()->SetBinLabel(2,"pt"); | |
614 | fhCutCorrelation->GetYaxis()->SetBinLabel(3,"px"); | |
615 | fhCutCorrelation->GetYaxis()->SetBinLabel(4,"py"); | |
616 | fhCutCorrelation->GetYaxis()->SetBinLabel(5,"pz"); | |
617 | fhCutCorrelation->GetYaxis()->SetBinLabel(6,"eta"); | |
618 | fhCutCorrelation->GetYaxis()->SetBinLabel(7,"y"); | |
619 | fhCutCorrelation->GetYaxis()->SetBinLabel(8,"phi"); | |
620 | fhCutCorrelation->GetYaxis()->SetBinLabel(9,"charge"); | |
621 | fhCutCorrelation->GetYaxis()->SetBinLabel(10,"is charged"); | |
622 | ||
623 | ||
624 | // book QA histograms | |
625 | Char_t str[256]; | |
626 | for (Int_t i=0; i<kNStepQA; i++) { | |
627 | if (i==0) sprintf(str," "); | |
628 | else sprintf(str,"_cut"); | |
629 | ||
630 | fhQA[kCutP][i] = new TH1F(Form("%s_momentum%s",GetName(),str), "",fhNBinsMomentum,fhBinLimMomentum); | |
631 | fhQA[kCutPt][i] = new TH1F(Form("%s_transverse_momentum%s",GetName(),str),"",fhNBinsPt,fhBinLimPt); | |
632 | fhQA[kCutPx][i] = new TH1F(Form("%s_px%s",GetName(),str), "",fhNBinsPx,fhBinLimPx); | |
633 | fhQA[kCutPy][i] = new TH1F(Form("%s_py%s",GetName(),str), "",fhNBinsPy,fhBinLimPy); | |
634 | fhQA[kCutPz][i] = new TH1F(Form("%s_pz%s",GetName(),str), "",fhNBinsPz,fhBinLimPz); | |
635 | fhQA[kCutRapidity][i]=new TH1F(Form("%s_rapidity%s",GetName(),str), "",fhNBinsRapidity,fhBinLimRapidity); | |
636 | fhQA[kCutEta][i] = new TH1F(Form("%s_eta%s",GetName(),str), "",fhNBinsEta,fhBinLimEta); | |
637 | fhQA[kCutPhi][i] = new TH1F(Form("%s_phi%s",GetName(),str), "",fhNBinsPhi,fhBinLimPhi); | |
638 | fhQA[kCutCharge][i] = new TH1F(Form("%s_charge%s",GetName(),str), "",fhNBinsCharge,fhBinLimCharge); | |
639 | ||
640 | fhQA[kCutP][i] ->SetXTitle("momentum p (GeV/c)"); | |
641 | fhQA[kCutPt][i] ->SetXTitle("p_{T} (GeV/c)"); | |
642 | fhQA[kCutPx][i] ->SetXTitle("p_{x} (GeV/c)"); | |
643 | fhQA[kCutPy][i] ->SetXTitle("p_{y} (GeV/c)"); | |
644 | fhQA[kCutPz][i] ->SetXTitle("p_{z} (GeV/c)"); | |
645 | fhQA[kCutRapidity][i]->SetXTitle("rapidity y"); | |
646 | fhQA[kCutEta][i] ->SetXTitle("pseudo rapidity #eta"); | |
647 | fhQA[kCutPhi][i] ->SetXTitle("azimuth #phi (rad)"); | |
648 | fhQA[kCutCharge][i] ->SetXTitle("charge"); | |
649 | } | |
650 | ||
651 | for(Int_t i=0; i<kNHist; i++) fhQA[i][1]->SetLineColor(color); | |
652 | ||
653 | } | |
654 | //__________________________________________________________________________________ | |
655 | void AliCFTrackKineCuts::FillHistograms(TObject* obj, Bool_t b) | |
656 | { | |
657 | // | |
658 | // fill the QA histograms | |
659 | // | |
660 | if(!fIsQAOn) return; | |
661 | ||
662 | // cast TObject into VParticle | |
663 | AliVParticle* particle = dynamic_cast<AliVParticle *>(obj); | |
664 | if ( !particle ) return; | |
665 | ||
666 | // index = 0: fill histograms before cuts | |
667 | // index = 1: fill histograms after cuts | |
668 | Int_t index = -1; | |
669 | index = ((b) ? 1 : 0); | |
670 | ||
671 | fhQA[kCutP][index]->Fill(particle->P()); | |
672 | fhQA[kCutPt][index]->Fill(particle->Pt()); | |
673 | fhQA[kCutPx][index]->Fill(particle->Px()); | |
674 | fhQA[kCutPy][index]->Fill(particle->Py()); | |
675 | fhQA[kCutPz][index]->Fill(particle->Pz()); | |
676 | fhQA[kCutRapidity][index]->Fill(particle->Y()); | |
677 | fhQA[kCutEta][index]->Fill(particle->Eta()); | |
678 | fhQA[kCutPhi][index]->Fill(particle->Phi()); | |
679 | fhQA[kCutCharge][index]->Fill((float)particle->Charge()); | |
680 | ||
681 | // fill cut statistics and cut correlation histograms with information from the bitmap | |
682 | if (b) return; | |
683 | ||
684 | if (!obj) return; | |
685 | TBits* bitmap = SelectionBitMap(obj); | |
686 | ||
687 | // Number of single cuts in this class | |
688 | UInt_t ncuts = bitmap->GetNbits(); | |
689 | for(UInt_t bit=0; bit<ncuts;bit++) { | |
690 | if (!bitmap->TestBitNumber(bit)) { | |
691 | fhCutStatistics->Fill(bit+1); | |
692 | for (UInt_t bit2=bit; bit2<ncuts;bit2++) { | |
693 | if (!bitmap->TestBitNumber(bit2)) | |
694 | fhCutCorrelation->Fill(bit+1,bit2+1); | |
695 | } | |
696 | } | |
697 | } | |
698 | } | |
699 | //__________________________________________________________________________________ | |
700 | void AliCFTrackKineCuts::SaveHistograms(const Char_t* dir) { | |
701 | // | |
702 | // saves the histograms in a directory (dir) | |
703 | // | |
704 | if(!fIsQAOn) return; | |
705 | ||
706 | if (!dir) | |
707 | dir = GetName(); | |
708 | ||
709 | gDirectory->mkdir(dir); | |
710 | gDirectory->cd(dir); | |
711 | ||
712 | gDirectory->mkdir("before_cuts"); | |
713 | gDirectory->mkdir("after_cuts"); | |
714 | ||
715 | fhCutStatistics->Write(); | |
716 | fhCutCorrelation->Write(); | |
717 | ||
718 | for (Int_t j=0; j<kNStepQA; j++) { | |
719 | if (j==0) | |
720 | gDirectory->cd("before_cuts"); | |
721 | else | |
722 | gDirectory->cd("after_cuts"); | |
723 | ||
724 | for(Int_t i=0; i<kNHist; i++) fhQA[i][j]->Write(); | |
725 | ||
726 | gDirectory->cd("../"); | |
727 | } | |
728 | ||
729 | gDirectory->cd("../"); | |
730 | } | |
731 | //__________________________________________________________________________________ | |
732 | void AliCFTrackKineCuts::DrawHistograms(Bool_t drawLogScale) | |
733 | { | |
734 | // | |
735 | // draws some histograms | |
736 | // | |
737 | if(!fIsQAOn) return; | |
738 | ||
739 | // pad margins | |
740 | Float_t right = 0.03; | |
741 | Float_t left = 0.175; | |
742 | Float_t top = 0.03; | |
743 | Float_t bottom = 0.175; | |
744 | ||
745 | TCanvas* canvas1 = new TCanvas("Track_QA_Kinematics_1", "Track QA Kinematics 1", 800, 500); | |
746 | canvas1->Divide(2, 1); | |
747 | ||
748 | canvas1->cd(1); | |
749 | fhCutStatistics->SetStats(kFALSE); | |
750 | fhCutStatistics->LabelsOption("v"); | |
751 | gPad->SetLeftMargin(left); | |
752 | gPad->SetBottomMargin(0.25); | |
753 | gPad->SetRightMargin(right); | |
754 | gPad->SetTopMargin(0.1); | |
755 | fhCutStatistics->Draw(); | |
756 | ||
757 | canvas1->cd(2); | |
758 | fhCutCorrelation->SetStats(kFALSE); | |
759 | fhCutCorrelation->LabelsOption("v"); | |
760 | gPad->SetLeftMargin(0.30); | |
761 | gPad->SetRightMargin(bottom); | |
762 | gPad->SetTopMargin(0.1); | |
763 | gPad->SetBottomMargin(0.25); | |
764 | fhCutCorrelation->Draw("COLZ"); | |
765 | ||
766 | canvas1->SaveAs(Form("%s.eps", canvas1->GetName())); | |
767 | canvas1->SaveAs(Form("%s.ps", canvas1->GetName())); | |
768 | ||
769 | // ----- | |
770 | ||
771 | TCanvas* canvas2 = new TCanvas("Track_QA_Kinematics_2", "Track QA Kinematics 2", 1600, 800); | |
772 | canvas2->Divide(4, 2); | |
773 | ||
774 | canvas2->cd(1); | |
775 | fhQA[kCutP][0]->SetStats(kFALSE); | |
776 | if(drawLogScale) gPad->SetLogy(); | |
777 | gPad->SetRightMargin(right); | |
778 | gPad->SetLeftMargin(left); | |
779 | gPad->SetTopMargin(top); | |
780 | gPad->SetBottomMargin(bottom); | |
781 | fhQA[kCutP][0]->Draw(); | |
782 | fhQA[kCutP][1]->Draw("same"); | |
783 | ||
784 | canvas2->cd(2); | |
785 | fhQA[kCutPt][0]->SetStats(kFALSE); | |
786 | if(drawLogScale) gPad->SetLogy(); | |
787 | gPad->SetRightMargin(right); | |
788 | gPad->SetLeftMargin(left); | |
789 | gPad->SetTopMargin(top); | |
790 | gPad->SetBottomMargin(bottom); | |
791 | fhQA[kCutPt][0]->Draw(); | |
792 | fhQA[kCutPt][1]->Draw("same"); | |
793 | ||
794 | canvas2->cd(3); | |
795 | fhQA[kCutRapidity][0]->SetStats(kFALSE); | |
796 | if(drawLogScale) gPad->SetLogy(); | |
797 | gPad->SetRightMargin(right); | |
798 | gPad->SetLeftMargin(left); | |
799 | gPad->SetTopMargin(top); | |
800 | gPad->SetBottomMargin(bottom); | |
801 | fhQA[kCutRapidity][0]->Draw(); | |
802 | fhQA[kCutRapidity][1]->Draw("same"); | |
803 | ||
804 | canvas2->cd(4); | |
805 | fhQA[kCutEta][0]->SetStats(kFALSE); | |
806 | if(drawLogScale) gPad->SetLogy(); | |
807 | gPad->SetRightMargin(right); | |
808 | gPad->SetLeftMargin(left); | |
809 | gPad->SetTopMargin(top); | |
810 | gPad->SetBottomMargin(bottom); | |
811 | fhQA[kCutEta][0]->Draw(); | |
812 | fhQA[kCutEta][1]->Draw("same"); | |
813 | ||
814 | canvas2->cd(5); | |
815 | fhQA[kCutPx][0]->SetStats(kFALSE); | |
816 | if(drawLogScale) gPad->SetLogy(); | |
817 | gPad->SetRightMargin(right); | |
818 | gPad->SetLeftMargin(left); | |
819 | gPad->SetTopMargin(top); | |
820 | gPad->SetBottomMargin(bottom); | |
821 | fhQA[kCutPx][0]->Draw(); | |
822 | fhQA[kCutPx][1]->Draw("same"); | |
823 | ||
824 | canvas2->cd(6); | |
825 | fhQA[kCutPy][0]->SetStats(kFALSE); | |
826 | if(drawLogScale) gPad->SetLogy(); | |
827 | gPad->SetRightMargin(right); | |
828 | gPad->SetLeftMargin(left); | |
829 | gPad->SetTopMargin(top); | |
830 | gPad->SetBottomMargin(bottom); | |
831 | fhQA[kCutPy][0]->Draw(); | |
832 | fhQA[kCutPy][1]->Draw("same"); | |
833 | ||
834 | canvas2->cd(7); | |
835 | fhQA[kCutPz][0]->SetStats(kFALSE); | |
836 | if(drawLogScale) gPad->SetLogy(); | |
837 | gPad->SetRightMargin(right); | |
838 | gPad->SetLeftMargin(left); | |
839 | gPad->SetTopMargin(top); | |
840 | gPad->SetBottomMargin(bottom); | |
841 | fhQA[kCutPz][0]->Draw(); | |
842 | fhQA[kCutPz][1]->Draw("same"); | |
843 | ||
844 | canvas2->SaveAs(Form("%s.eps", canvas2->GetName())); | |
845 | canvas2->SaveAs(Form("%s.ps", canvas2->GetName())); | |
846 | ||
847 | // ----- | |
848 | ||
849 | TCanvas* canvas3 = new TCanvas("Track_QA_Kinematics_3", "Track QA Kinematics 3", 800, 400); | |
850 | canvas3->Divide(2, 1); | |
851 | ||
852 | canvas3->cd(1); | |
853 | gPad->SetRightMargin(right); | |
854 | gPad->SetLeftMargin(left); | |
855 | gPad->SetTopMargin(top); | |
856 | gPad->SetBottomMargin(bottom); | |
857 | fhQA[kCutPhi][0]->SetStats(kFALSE); | |
858 | fhQA[kCutPhi][0]->Draw(); | |
859 | fhQA[kCutPhi][1]->Draw("same"); | |
860 | ||
861 | canvas3->cd(2); | |
862 | gPad->SetRightMargin(right); | |
863 | gPad->SetLeftMargin(left); | |
864 | gPad->SetTopMargin(top); | |
865 | gPad->SetBottomMargin(bottom); | |
866 | fhQA[kCutCharge][0]->SetStats(kFALSE); | |
867 | fhQA[kCutCharge][0]->Draw(); | |
868 | fhQA[kCutCharge][1]->Draw("same"); | |
869 | ||
870 | canvas3->SaveAs(Form("%s.eps", canvas3->GetName())); | |
871 | canvas3->SaveAs(Form("%s.ps", canvas3->GetName())); | |
872 | } | |
873 | //__________________________________________________________________________________ | |
874 | void AliCFTrackKineCuts::AddQAHistograms(TList *qaList) const { | |
875 | // | |
876 | // saves the histograms in a TList | |
877 | // | |
878 | if(!fIsQAOn) return; | |
879 | ||
880 | qaList->Add(fhCutStatistics); | |
881 | qaList->Add(fhCutCorrelation); | |
882 | ||
883 | for (Int_t j=0; j<kNStepQA; j++) { | |
884 | for(Int_t i=0; i<kNHist; i++) | |
885 | qaList->Add(fhQA[i][j]); | |
886 | } | |
887 | } |