]> git.uio.no Git - u/mrichter/AliRoot.git/blob - CORRFW/AliCFTrackIsPrimaryCuts.cxx
Macro to produce SPD vertices with Z and 3D vertexers
[u/mrichter/AliRoot.git] / CORRFW / AliCFTrackIsPrimaryCuts.cxx
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 AliCFTrackIsPrimaryCut is designed to select reconstructed tracks
17 // with a small impact parameter and tracks which are (not) daughters of kink
18 // decays and to provide corresponding QA histograms.
19 // This class inherits from the Analysis' Framework abstract base class
20 // AliAnalysisCuts and is a part of the Correction Framework.
21 // This class acts on single, reconstructed tracks, it is applicable on
22 // ESD and AOD data.
23 // It mainly consists of a IsSelected function that returns a boolean.
24 // This function checks whether the considered track passes a set of cuts:
25 // - distance to main vertex in units of sigma (resolution)
26 // - require that the dca calculation doesn't fail
27 // - accept or not accept daughter tracks of kink decays
28 //
29 // The cut values for these cuts are set with the corresponding set functions.
30 // All cut classes provided by the correction framework are supposed to be
31 // added in the Analysis Framwork's class AliAnalysisFilter and applied by
32 // the filter via a loop.
33 //
34 // author: I. Kraus (Ingrid.Kraus@cern.ch)
35 // idea taken form
36 // AliESDtrackCuts writte by Jan Fiete Grosse-Oetringhaus and
37 // AliRsnDaughterCut class written by A. Pulvirenti.
38
39 #include <TCanvas.h>
40 #include <TDirectory.h>
41 #include <TH2.h>
42 #include <TBits.h>
43
44 #include <AliESDtrack.h>
45 #include <AliLog.h>
46 #include "AliCFTrackIsPrimaryCuts.h"
47
48 ClassImp(AliCFTrackIsPrimaryCuts)
49
50 //__________________________________________________________________________________
51 AliCFTrackIsPrimaryCuts::AliCFTrackIsPrimaryCuts() :
52   AliCFCutBase(),
53   fNSigmaToVertex(0),
54   fNSigmaToVertexMax(0),
55   fRequireSigmaToVertex(0),
56   fAcceptKinkDaughters(0),
57   fhCutStatistics(0),
58   fhCutCorrelation(0),
59   fBitmap(0x0),
60   fhNBinsNSigma(0),
61   fhNBinsRequireSigma(0),
62   fhNBinsAcceptKink(0),
63   fhNBinsDcaXY(0),
64   fhNBinsDcaZ(0),
65   fhNBinsDcaXYnorm(0),
66   fhNBinsDcaZnorm(0),
67   fhBinLimNSigma(0x0),
68   fhBinLimRequireSigma(0x0),
69   fhBinLimAcceptKink(0x0),
70   fhBinLimDcaXY(0x0),
71   fhBinLimDcaZ(0x0),
72   fhBinLimDcaXYnorm(0x0),
73   fhBinLimDcaZnorm(0x0)
74 {
75   //
76   // Default constructor
77   //
78   Initialise();
79 }
80 //__________________________________________________________________________________
81 AliCFTrackIsPrimaryCuts::AliCFTrackIsPrimaryCuts(Char_t* name, Char_t* title) :
82   AliCFCutBase(name,title),
83   fNSigmaToVertex(0),
84   fNSigmaToVertexMax(0),
85   fRequireSigmaToVertex(0),
86   fAcceptKinkDaughters(0),
87   fhCutStatistics(0),
88   fhCutCorrelation(0),
89   fBitmap(0x0),
90   fhNBinsNSigma(0),
91   fhNBinsRequireSigma(0),
92   fhNBinsAcceptKink(0),
93   fhNBinsDcaXY(0),
94   fhNBinsDcaZ(0),
95   fhNBinsDcaXYnorm(0),
96   fhNBinsDcaZnorm(0),
97   fhBinLimNSigma(0x0),
98   fhBinLimRequireSigma(0x0),
99   fhBinLimAcceptKink(0x0),
100   fhBinLimDcaXY(0x0),
101   fhBinLimDcaZ(0x0),
102   fhBinLimDcaXYnorm(0x0),
103   fhBinLimDcaZnorm(0x0)
104 {
105   //
106   // Constructor
107   //
108   Initialise();
109 }
110 //__________________________________________________________________________________
111 AliCFTrackIsPrimaryCuts::AliCFTrackIsPrimaryCuts(const AliCFTrackIsPrimaryCuts& c) :
112   AliCFCutBase(c),
113   fNSigmaToVertex(c.fNSigmaToVertex),
114   fNSigmaToVertexMax(c.fNSigmaToVertexMax),
115   fRequireSigmaToVertex(c.fRequireSigmaToVertex),
116   fAcceptKinkDaughters(c.fAcceptKinkDaughters),
117   fhCutStatistics(c.fhCutStatistics),
118   fhCutCorrelation(c.fhCutCorrelation),
119   fBitmap(c.fBitmap),
120   fhNBinsNSigma(c.fhNBinsNSigma),
121   fhNBinsRequireSigma(c.fhNBinsRequireSigma),
122   fhNBinsAcceptKink(c.fhNBinsAcceptKink),
123   fhNBinsDcaXY(c.fhNBinsDcaXY),
124   fhNBinsDcaZ(c.fhNBinsDcaZ),
125   fhNBinsDcaXYnorm(c.fhNBinsDcaXYnorm),
126   fhNBinsDcaZnorm(c.fhNBinsDcaZnorm),
127   fhBinLimNSigma(c.fhBinLimNSigma),
128   fhBinLimRequireSigma(c.fhBinLimRequireSigma),
129   fhBinLimAcceptKink(c.fhBinLimAcceptKink),
130   fhBinLimDcaXY(c.fhBinLimDcaXY),
131   fhBinLimDcaZ(c.fhBinLimDcaZ),
132   fhBinLimDcaXYnorm(c.fhBinLimDcaXYnorm),
133   fhBinLimDcaZnorm(c.fhBinLimDcaZnorm)
134 {
135   //
136   // copy constructor
137   //
138   ((AliCFTrackIsPrimaryCuts &) c).Copy(*this);
139 }
140 //__________________________________________________________________________________
141 AliCFTrackIsPrimaryCuts& AliCFTrackIsPrimaryCuts::operator=(const AliCFTrackIsPrimaryCuts& c)
142 {
143   //
144   // Assignment operator
145   //
146   if (this != &c) {
147     AliCFCutBase::operator=(c) ;
148     fNSigmaToVertex = c.fNSigmaToVertex ;
149     fNSigmaToVertexMax = c.fNSigmaToVertexMax ;
150     fRequireSigmaToVertex = c.fRequireSigmaToVertex ;
151     fAcceptKinkDaughters = c.fAcceptKinkDaughters ;
152     fhCutStatistics = c.fhCutStatistics ;
153     fhCutCorrelation = c.fhCutCorrelation ;
154     fBitmap =  c.fBitmap;
155     fhNBinsNSigma = c.fhNBinsNSigma;
156     fhNBinsRequireSigma = c.fhNBinsRequireSigma;
157     fhNBinsAcceptKink = c.fhNBinsAcceptKink;
158     fhNBinsDcaXY = c.fhNBinsDcaXY;
159     fhNBinsDcaZ = c.fhNBinsDcaZ;
160     fhNBinsDcaXYnorm = c.fhNBinsDcaXYnorm;
161     fhNBinsDcaZnorm = c.fhNBinsDcaZnorm;
162     fhBinLimNSigma = c.fhBinLimNSigma;
163     fhBinLimRequireSigma = c.fhBinLimRequireSigma;
164     fhBinLimAcceptKink = c.fhBinLimAcceptKink;
165     fhBinLimDcaXY = c.fhBinLimDcaXY;
166     fhBinLimDcaZ = c.fhBinLimDcaZ;
167     fhBinLimDcaXYnorm = c.fhBinLimDcaXYnorm;
168     fhBinLimDcaZnorm = c.fhBinLimDcaZnorm;
169
170     for (Int_t i=0; i<c.kNHist; i++){
171       for (Int_t j=0; j<c.kNStepQA; j++){
172         if(c.fhQA[i][j]) fhQA[i][j] = (TH1F*)c.fhQA[i][j]->Clone();
173         if(c.fhDcaXYvsDcaZ[j]) fhDcaXYvsDcaZ[j] = (TH2F*)c.fhDcaXYvsDcaZ[j]->Clone();
174         if(c.fhDcaXYvsDcaZnorm[j]) fhDcaXYvsDcaZnorm[j] = (TH2F*)c.fhDcaXYvsDcaZnorm[j]->Clone();
175       }
176     }
177     ((AliCFTrackIsPrimaryCuts &) c).Copy(*this);
178  }
179   return *this;
180 }
181 //__________________________________________________________________________________
182 AliCFTrackIsPrimaryCuts::~AliCFTrackIsPrimaryCuts()
183 {
184   //
185   // destructor
186   //
187   if (fhCutStatistics)                  delete fhCutStatistics;
188   if (fhCutCorrelation)                 delete fhCutCorrelation;
189
190   for (Int_t j=0; j<kNStepQA; j++){
191     if(fhDcaXYvsDcaZ[j])        delete fhDcaXYvsDcaZ[j];
192     if(fhDcaXYvsDcaZnorm[j])    delete fhDcaXYvsDcaZnorm[j];
193     for (Int_t i=0; i<kNHist; i++)
194       if(fhQA[i][j])            delete fhQA[i][j];
195   }
196   if(fBitmap)                   delete fBitmap;
197   if(fhBinLimNSigma)            delete fhBinLimNSigma;
198   if(fhBinLimRequireSigma)      delete fhBinLimRequireSigma;
199   if(fhBinLimAcceptKink)        delete fhBinLimAcceptKink;
200   if(fhBinLimDcaXY)             delete fhBinLimDcaXY;
201   if(fhBinLimDcaZ)              delete fhBinLimDcaZ;
202   if(fhBinLimDcaXYnorm)         delete fhBinLimDcaXYnorm;
203   if(fhBinLimDcaZnorm)          delete fhBinLimDcaZnorm;
204 }
205 //__________________________________________________________________________________
206 void AliCFTrackIsPrimaryCuts::Initialise()
207 {
208   //
209   // sets everything to zero
210   //
211   fNSigmaToVertex = 0;
212   fNSigmaToVertexMax = 0;
213   fRequireSigmaToVertex = 0;
214   fAcceptKinkDaughters = 0;
215
216   SetMaxNSigmaToVertex();
217   SetRequireSigmaToVertex();
218
219   for (Int_t j=0; j<kNStepQA; j++)  {
220     fhDcaXYvsDcaZ[j] = 0x0;
221     fhDcaXYvsDcaZnorm[j] = 0x0;
222     for (Int_t i=0; i<kNHist; i++)
223       fhQA[i][j] = 0x0;
224   }
225   fhCutStatistics = 0;
226   fhCutCorrelation = 0;
227   fBitmap=new TBits(0);
228
229   //set default bining for QA histograms
230   SetHistogramBins(kCutNSigmaToVertex,500,0.,50.);
231   SetHistogramBins(kCutRequireSigmaToVertex,5,-0.75,1.75);
232   SetHistogramBins(kCutAcceptKinkDaughters,5,-0.75,1.75);
233   SetHistogramBins(kDcaXY,500,-10.,10.);
234   SetHistogramBins(kDcaZ,500,-10.,10.);
235   SetHistogramBins(kDcaXYnorm,500,-10.,10.);
236   SetHistogramBins(kDcaZnorm,500,-10.,10.);
237 }
238 //__________________________________________________________________________________
239 void AliCFTrackIsPrimaryCuts::Copy(TObject &c) const
240 {
241   //
242   // Copy function
243   //
244   AliCFTrackIsPrimaryCuts& target = (AliCFTrackIsPrimaryCuts &) c;
245
246   target.Initialise();
247
248   if (fhCutStatistics)  target.fhCutStatistics = (TH1F*) fhCutStatistics->Clone();
249   if (fhCutCorrelation) target.fhCutCorrelation = (TH2F*) fhCutCorrelation->Clone();
250
251   for (Int_t j=0; j<kNStepQA; j++){
252     if(fhDcaXYvsDcaZ[j]) target.fhDcaXYvsDcaZ[j] = (TH2F*)fhDcaXYvsDcaZ[j]->Clone();
253     if(fhDcaXYvsDcaZnorm[j]) target.fhDcaXYvsDcaZnorm[j] = (TH2F*)fhDcaXYvsDcaZnorm[j]->Clone();
254     for (Int_t i=0; i<kNHist; i++)
255       if(fhQA[i][j]) target.fhQA[i][j] = (TH1F*)fhQA[i][j]->Clone();
256   }
257   TNamed::Copy(c);
258 }
259 //____________________________________________________________________
260 void AliCFTrackIsPrimaryCuts::GetSigmaToVertex(AliESDtrack* esdTrack)
261 {
262   //
263   // Calculates the number of sigma to the vertex.
264   //
265   Float_t b[2];
266   Float_t bRes[2];
267   Float_t bCov[3];
268   esdTrack->GetImpactParameters(b,bCov);
269   if (bCov[0]<=0 || bCov[2]<=0) {
270     AliDebug(1, "Estimated b resolution lower or equal zero!");
271     bCov[0]=0; bCov[2]=0;
272   }
273   bRes[0] = TMath::Sqrt(bCov[0]);
274   bRes[1] = TMath::Sqrt(bCov[2]);
275
276   // -----------------------------------
277   // How to get to a n-sigma cut?
278   //
279   // The accumulated statistics from 0 to d is
280   //
281   // ->  Erf(d/Sqrt(2)) for a 1-dim gauss (d = n_sigma)
282   // ->  1 - Exp(-d**2) for a 2-dim gauss (d*d = dx*dx + dy*dy != n_sigma)
283   //
284   // It means that for a 2-dim gauss: n_sigma(d) = Sqrt(2)*ErfInv(1 - Exp((-x**2)/2)
285   // Can this be expressed in a different way?
286
287   if (bRes[0] == 0 || bRes[1] ==0) {
288     fNSigmaToVertex = -1;
289     return;
290   }
291
292   Float_t d = TMath::Sqrt(TMath::Power(b[0]/bRes[0],2) + TMath::Power(b[1]/bRes[1],2));
293
294   // stupid rounding problem screws up everything:
295   // if d is too big, TMath::Exp(...) gets 0, and TMath::ErfInverse(1) that should be infinite, gets 0 :(
296   if (TMath::Exp(-d * d / 2) < 1e-10) {
297     fNSigmaToVertex = 1000;
298     return;
299   }
300
301   d = TMath::ErfInverse(1 - TMath::Exp(-d * d / 2)) * TMath::Sqrt(2);
302   fNSigmaToVertex = d;
303   return;
304 }
305 //__________________________________________________________________________________
306 void AliCFTrackIsPrimaryCuts::SelectionBitMap(TObject* obj)
307 {
308   //
309   // test if the track passes the single cuts
310   // and store the information in a bitmap
311   //
312
313   // bitmap stores the decision of each single cut
314   for(Int_t i=0; i<kNCuts; i++)fBitmap->SetBitNumber(i,kFALSE);
315
316   // check TObject and cast into ESDtrack
317   if (!obj) return;
318   if (!obj->InheritsFrom("AliESDtrack")) AliError("object must derived from AliESDtrack !");
319   AliESDtrack* esdTrack = dynamic_cast<AliESDtrack *>(obj);
320   if ( !esdTrack ) return;
321
322   // get the track to vertex parameter
323   GetSigmaToVertex(esdTrack);
324
325   // fill the bitmap
326   Int_t iCutBit = 0;
327   if (fNSigmaToVertex <= fNSigmaToVertexMax)
328         fBitmap->SetBitNumber(iCutBit,kTRUE);
329   iCutBit++;
330   if (!fRequireSigmaToVertex || (fNSigmaToVertex>=0 && fRequireSigmaToVertex))
331         fBitmap->SetBitNumber(iCutBit,kTRUE);
332   iCutBit++;
333   if (fAcceptKinkDaughters || (!fAcceptKinkDaughters && esdTrack->GetKinkIndex(0)<=0))
334         fBitmap->SetBitNumber(iCutBit,kTRUE);
335
336   return;
337 }
338 //__________________________________________________________________________________
339 Bool_t AliCFTrackIsPrimaryCuts::IsSelected(TObject* obj) {
340   //
341   // loops over decisions of single cuts and returns if the track is accepted
342   //
343   SelectionBitMap(obj);
344
345   if (fIsQAOn) FillHistograms(obj,0);
346   Bool_t isSelected = kTRUE;
347
348   for (UInt_t icut=0; icut<fBitmap->GetNbits();icut++) {
349     if(!fBitmap->TestBitNumber(icut)) {
350         isSelected = kFALSE;
351         break;
352     }
353   }
354   if (!isSelected) return kFALSE ;
355   if (fIsQAOn) FillHistograms(obj,1);
356   return kTRUE;
357 }
358 //__________________________________________________________________________________
359 void AliCFTrackIsPrimaryCuts::SetHistogramBins(Int_t index, Int_t nbins, Double_t *bins)
360 {
361   //
362   // variable bin size
363   //
364   if(!fIsQAOn) return;
365
366   switch(index){
367   case kCutNSigmaToVertex:
368     fhNBinsNSigma=nbins+1;
369     fhBinLimNSigma=new Double_t[nbins+1];
370     for(Int_t i=0;i<nbins+1;i++)fhBinLimNSigma[i]=bins[i];
371     break;
372
373   case kCutRequireSigmaToVertex:
374     fhNBinsRequireSigma=nbins+1;
375     fhBinLimRequireSigma=new Double_t[nbins+1];
376     for(Int_t i=0;i<nbins+1;i++)fhBinLimRequireSigma[i]=bins[i];
377     break;
378
379   case kCutAcceptKinkDaughters:
380     fhNBinsAcceptKink=nbins+1;
381     fhBinLimAcceptKink=new Double_t[nbins+1];
382     for(Int_t i=0;i<nbins+1;i++)fhBinLimAcceptKink[i]=bins[i];
383     break;
384
385   case kDcaXY:
386     fhNBinsDcaXY=nbins+1;
387     fhBinLimDcaXY=new Double_t[nbins+1];
388     for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaXY[i]=bins[i];
389     break;
390
391   case kDcaZ:
392     fhNBinsDcaZ=nbins+1;
393     fhBinLimDcaZ=new Double_t[nbins+1];
394     for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaZ[i]=bins[i];
395     break;
396
397   case kDcaXYnorm:
398     fhNBinsDcaXYnorm=nbins+1;
399     fhBinLimDcaXYnorm=new Double_t[nbins+1];
400     for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaXYnorm[i]=bins[i];
401     break;
402
403   case kDcaZnorm:
404     fhNBinsDcaZnorm=nbins+1;
405     fhBinLimDcaZnorm=new Double_t[nbins+1];
406     for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaZnorm[i]=bins[i];
407     break;
408   }
409 }
410 //__________________________________________________________________________________
411 void AliCFTrackIsPrimaryCuts::SetHistogramBins(Int_t index, Int_t nbins, Double_t xmin, Double_t xmax)
412 {
413   //
414   // fixed bin size
415   //
416   switch(index){
417   case kCutNSigmaToVertex:
418     fhNBinsNSigma=nbins+1;
419     fhBinLimNSigma=new Double_t[nbins+1];
420     for(Int_t i=0;i<nbins+1;i++)fhBinLimNSigma[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
421     break;
422
423   case kCutRequireSigmaToVertex:
424     fhNBinsRequireSigma=nbins+1;
425     fhBinLimRequireSigma=new Double_t[nbins+1];
426     for(Int_t i=0;i<nbins+1;i++)fhBinLimRequireSigma[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
427     break;
428
429   case kCutAcceptKinkDaughters:
430     fhNBinsAcceptKink=nbins+1;
431     fhBinLimAcceptKink=new Double_t[nbins+1];
432     for(Int_t i=0;i<nbins+1;i++)fhBinLimAcceptKink[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
433     break;
434
435   case kDcaXY:
436     fhNBinsDcaXY=nbins+1;
437     fhBinLimDcaXY=new Double_t[nbins+1];
438     for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaXY[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
439     break;
440
441   case kDcaZ:
442     fhNBinsDcaZ=nbins+1;
443     fhBinLimDcaZ=new Double_t[nbins+1];
444     for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaZ[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
445     break;
446
447   case kDcaXYnorm:
448     fhNBinsDcaXYnorm=nbins+1;
449     fhBinLimDcaXYnorm=new Double_t[nbins+1];
450     for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaXYnorm[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
451     break;
452
453   case kDcaZnorm:
454     fhNBinsDcaZnorm=nbins+1;
455     fhBinLimDcaZnorm=new Double_t[nbins+1];
456     for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaZnorm[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
457     break;
458   }
459 }
460 //__________________________________________________________________________________
461  void AliCFTrackIsPrimaryCuts::DefineHistograms() {
462   //
463   // histograms for cut variables, cut statistics and cut correlations
464   //
465   Int_t color = 2;
466
467   // book cut statistics and cut correlation histograms
468   fhCutStatistics = new TH1F(Form("%s_cut_statistics",GetName()), Form("%s cut statistics",GetName()), kNCuts,0.5,kNCuts+0.5);
469   fhCutStatistics->SetLineWidth(2);
470   fhCutStatistics->GetXaxis()->SetBinLabel(1,"n dca");
471   fhCutStatistics->GetXaxis()->SetBinLabel(2,"require dca");
472   fhCutStatistics->GetXaxis()->SetBinLabel(3,"kink daughter");
473
474   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);
475   fhCutCorrelation->SetLineWidth(2);
476   fhCutCorrelation->GetXaxis()->SetBinLabel(1,"n dca");
477   fhCutCorrelation->GetXaxis()->SetBinLabel(2,"require dca");
478   fhCutCorrelation->GetXaxis()->SetBinLabel(3,"kink daughter");
479
480   fhCutCorrelation->GetYaxis()->SetBinLabel(1,"n dca");
481   fhCutCorrelation->GetYaxis()->SetBinLabel(2,"require dca");
482   fhCutCorrelation->GetYaxis()->SetBinLabel(3,"kink daughter");
483
484   // book QA histograms
485   Char_t str[256];
486   for (Int_t i=0; i<kNStepQA; i++) {
487     if (i==0) sprintf(str," ");
488     else sprintf(str,"_cut");
489
490     fhDcaXYvsDcaZ[i]            = new  TH2F(Form("%s_dcaXYvsDcaZ%s",GetName(),str),"",200,-10,10,200,-10,10);
491     fhDcaXYvsDcaZnorm[i]        = new  TH2F(Form("%s_dcaXYvsDcaZnorm%s",GetName(),str),"",200,-10,10,200,-10,10);
492     fhQA[kCutNSigmaToVertex][i] = new TH1F(Form("%s_nSigmaToVertex%s",GetName(),str),"",fhNBinsNSigma-1,fhBinLimNSigma);
493     fhQA[kCutRequireSigmaToVertex][i] = new TH1F(Form("%s_requireSigmaToVertex%s",GetName(),str),"",fhNBinsRequireSigma-1,fhBinLimRequireSigma);
494     fhQA[kCutAcceptKinkDaughters][i] = new TH1F(Form("%s_acceptKinkDaughters%s",GetName(),str),"",fhNBinsAcceptKink-1,fhBinLimAcceptKink);
495     fhQA[kDcaXY][i]             = new TH1F(Form("%s_dcaXY%s",GetName(),str),"",fhNBinsDcaXY-1,fhBinLimDcaXY);
496     fhQA[kDcaZ][i]              = new TH1F(Form("%s_dcaZ%s",GetName(),str),"",fhNBinsDcaZ-1,fhBinLimDcaZ);
497     fhQA[kDcaXYnorm][i]         = new TH1F(Form("%s_dcaXYnorm%s",GetName(),str),"",fhNBinsDcaXYnorm-1,fhBinLimDcaXYnorm);
498     fhQA[kDcaZnorm][i]          = new TH1F(Form("%s_dcaZnorm%s",GetName(),str),"",fhNBinsDcaZnorm-1,fhBinLimDcaZnorm);
499
500     fhDcaXYvsDcaZ[i]->SetXTitle("impact par. d_{z}");
501     fhDcaXYvsDcaZ[i]->SetYTitle("impact par. d_{xy}");
502     fhDcaXYvsDcaZnorm[i]->SetXTitle("norm. impact par. d_{z} / #sigma_{z}");
503     fhDcaXYvsDcaZnorm[i]->SetYTitle("norm. impact par. d_{xy} / #sigma_{xy}");
504
505     fhQA[kCutNSigmaToVertex][i]->SetXTitle("n #sigma to vertex");
506     fhQA[kCutRequireSigmaToVertex][i]->SetXTitle("require #sigma to vertex");
507     fhQA[kCutAcceptKinkDaughters][i]->SetXTitle("accept kink daughters");
508     fhQA[kDcaXY][i]->SetXTitle("impact par. d_{xy}");
509     fhQA[kDcaZ][i]->SetXTitle("impact par. d_{z}");
510     fhQA[kDcaXYnorm][i]->SetXTitle("norm. impact par. d_{xy} / #sigma_{xy}");
511     fhQA[kDcaZnorm][i]->SetXTitle("norm. impact par. d_{z} / #sigma_{z}");
512   }
513   for(Int_t i=0; i<kNHist; i++) fhQA[i][1]->SetLineColor(color);
514 }
515 //__________________________________________________________________________________
516 void AliCFTrackIsPrimaryCuts::FillHistograms(TObject* obj, Bool_t f)
517 {
518   //
519   // fill the QA histograms
520   //
521
522   // cast TObject into ESDtrack
523   if (!obj) return;
524   AliESDtrack* esdTrack = dynamic_cast<AliESDtrack *>(obj);
525   if ( !esdTrack ) return;
526
527   // f = 0: fill histograms before cuts
528   // f = 1: fill histograms after cuts
529
530         Float_t b[2];
531         Float_t bRes[2];
532         Float_t bCov[3];
533         esdTrack->GetImpactParameters(b,bCov);
534         if (bCov[0]<=0 || bCov[2]<=0) {
535            AliDebug(1, "Estimated b resolution lower or equal zero!");
536            bCov[0]=0; bCov[2]=0;
537         }
538         bRes[0] = TMath::Sqrt(bCov[0]);
539         bRes[1] = TMath::Sqrt(bCov[2]);
540
541         fhQA[kDcaZ][f]->Fill(b[1]);
542         fhQA[kDcaXY][f]->Fill(b[0]);
543         fhDcaXYvsDcaZ[f]->Fill(b[1],b[0]);
544
545         if (bRes[0]!=0 && bRes[1]!=0) {
546            fhQA[kDcaZnorm][f]->Fill(b[1]/bRes[1]);
547            fhQA[kDcaXYnorm][f]->Fill(b[0]/bRes[0]);
548            fhDcaXYvsDcaZnorm[f]->Fill(b[1]/bRes[1], b[0]/bRes[0]);
549         }
550
551   fhQA[kCutNSigmaToVertex][f]->Fill(fNSigmaToVertex);
552   if (fNSigmaToVertex<0 && fRequireSigmaToVertex) fhQA[kCutRequireSigmaToVertex][f]->Fill(0.);
553   if (!(fNSigmaToVertex<0 && fRequireSigmaToVertex)) fhQA[kCutRequireSigmaToVertex][f]->Fill(1.);
554
555
556   if (!fAcceptKinkDaughters && esdTrack->GetKinkIndex(0)>0) fhQA[kCutAcceptKinkDaughters][f]->Fill(0.);
557   if (!(!fAcceptKinkDaughters && esdTrack->GetKinkIndex(0)>0)) fhQA[kCutAcceptKinkDaughters][f]->Fill(0.);
558
559   // fill cut statistics and cut correlation histograms with information from the bitmap
560   if (f) return;
561
562   // Get the bitmap of the single cuts
563   if ( !obj ) return;
564   SelectionBitMap(obj);
565
566   // Number of single cuts in this class
567   UInt_t ncuts = fBitmap->GetNbits();
568   for(UInt_t bit=0; bit<ncuts;bit++) {
569     if (!fBitmap->TestBitNumber(bit)) {
570         fhCutStatistics->Fill(bit+1);
571         for (UInt_t bit2=bit; bit2<ncuts;bit2++) {
572           if (!fBitmap->TestBitNumber(bit2)) 
573             fhCutCorrelation->Fill(bit+1,bit2+1);
574         }
575     }
576   }
577 }
578 //__________________________________________________________________________________
579 void AliCFTrackIsPrimaryCuts::SaveHistograms(const Char_t* dir) {
580   //
581   // saves the histograms in a directory (dir)
582   //
583   if(!fIsQAOn) return;
584
585   if (!dir)
586     dir = GetName();
587
588   gDirectory->mkdir(dir);
589   gDirectory->cd(dir);
590
591   gDirectory->mkdir("before_cuts");
592   gDirectory->mkdir("after_cuts");
593
594   fhCutStatistics->Write();
595   fhCutCorrelation->Write();
596
597   for (Int_t j=0; j<kNStepQA; j++) {
598     if (j==0)
599       gDirectory->cd("before_cuts");
600     else
601       gDirectory->cd("after_cuts");
602
603     fhDcaXYvsDcaZ[j]    ->Write();
604     fhDcaXYvsDcaZnorm[j]->Write();
605
606     for(Int_t i=0; i<kNHist; i++) fhQA[i][j]->Write();
607
608     gDirectory->cd("../");
609   }
610   gDirectory->cd("../");
611 }
612 //__________________________________________________________________________________
613 void AliCFTrackIsPrimaryCuts::DrawHistograms()
614 {
615   //
616   // draws some histograms
617   //
618   if(!fIsQAOn) return;
619
620   // pad margins
621   Float_t right = 0.03;
622   Float_t left = 0.175;
623   Float_t top = 0.03;
624   Float_t bottom = 0.175;
625
626   TCanvas* canvas1 = new TCanvas("Track_QA_Primary_1", "Track QA Primary 1", 800, 500);
627   canvas1->Divide(2, 1);
628
629   canvas1->cd(1);
630   fhCutStatistics->SetStats(kFALSE);
631   fhCutStatistics->LabelsOption("v");
632   gPad->SetLeftMargin(left);
633   gPad->SetBottomMargin(0.25);
634   gPad->SetRightMargin(right);
635   gPad->SetTopMargin(0.1);
636   fhCutStatistics->Draw();
637
638   canvas1->cd(2);
639   fhCutCorrelation->SetStats(kFALSE);
640   fhCutCorrelation->LabelsOption("v");
641   gPad->SetLeftMargin(0.30);
642   gPad->SetRightMargin(bottom);
643   gPad->SetTopMargin(0.1);
644   gPad->SetBottomMargin(0.25);
645   fhCutCorrelation->Draw("COLZ");
646
647   canvas1->SaveAs(Form("%s.eps", canvas1->GetName()));
648   canvas1->SaveAs(Form("%s.ps", canvas1->GetName()));
649
650   // -----
651
652   TCanvas* canvas2 = new TCanvas("Track_QA_Primary_2", "Track QA Primary 2", 800, 800);
653   canvas2->Divide(2, 2);
654
655   canvas2->cd(1);
656   gPad->SetRightMargin(right);
657   gPad->SetLeftMargin(left);
658   gPad->SetTopMargin(top);
659   gPad->SetBottomMargin(bottom);
660   fhQA[kDcaXY][0]->SetStats(kFALSE);
661   fhQA[kDcaXY][0]->Draw();
662   fhQA[kDcaXY][1]->Draw("same");
663
664   canvas2->cd(2);
665   gPad->SetRightMargin(right);
666   gPad->SetLeftMargin(left);
667   gPad->SetTopMargin(top);
668   gPad->SetBottomMargin(bottom);
669   fhQA[kDcaZ][0]->SetStats(kFALSE);
670   fhQA[kDcaZ][0]->Draw();
671   fhQA[kDcaZ][1]->Draw("same");
672
673   canvas2->cd(3);
674 //   fhDXYvsDZ[0]->SetStats(kFALSE);
675   gPad->SetLogz();
676   gPad->SetLeftMargin(bottom);
677   gPad->SetTopMargin(0.1);
678   gPad->SetBottomMargin(bottom);
679   gPad->SetRightMargin(0.2);
680   fhDcaXYvsDcaZ[0]->Draw("COLZ");
681
682   canvas2->cd(4);
683 //   fhDXYvsDZ[1]->SetStats(kFALSE);
684   gPad->SetLogz();
685   gPad->SetLeftMargin(bottom);
686   gPad->SetTopMargin(0.1);
687   gPad->SetBottomMargin(bottom);
688   gPad->SetRightMargin(0.2);
689   fhDcaXYvsDcaZ[1]->Draw("COLZ");
690
691   canvas2->SaveAs(Form("%s.eps", canvas2->GetName()));
692   canvas2->SaveAs(Form("%s.ps", canvas2->GetName()));
693
694   // -----
695
696   TCanvas* canvas3 = new TCanvas("Track_QA_Primary_3", "Track QA Primary 3", 800, 800);
697   canvas3->Divide(2, 2);
698
699   canvas3->cd(1);
700   gPad->SetRightMargin(right);
701   gPad->SetLeftMargin(left);
702   gPad->SetTopMargin(top);
703   gPad->SetBottomMargin(bottom);
704   fhQA[kDcaXYnorm][0]->SetStats(kFALSE);
705   fhQA[kDcaXYnorm][0]->Draw();
706   fhQA[kDcaXYnorm][1]->Draw("same");
707
708   canvas3->cd(2);
709   gPad->SetRightMargin(right);
710   gPad->SetLeftMargin(left);
711   gPad->SetTopMargin(top);
712   gPad->SetBottomMargin(bottom);
713   fhQA[kDcaZnorm][0]->SetStats(kFALSE);
714   fhQA[kDcaZnorm][0]->Draw();
715   fhQA[kDcaZnorm][1]->Draw("same");
716
717   canvas3->cd(3);
718 //   fhDXYvsDZ[0]->SetStats(kFALSE);
719   gPad->SetLogz();
720   gPad->SetLeftMargin(bottom);
721   gPad->SetTopMargin(0.1);
722   gPad->SetBottomMargin(bottom);
723   gPad->SetRightMargin(0.2);
724   fhDcaXYvsDcaZnorm[0]->Draw("COLZ");
725
726   canvas3->cd(4);
727 //   fhDXYvsDZ[1]->SetStats(kFALSE);
728   gPad->SetLogz();
729   gPad->SetLeftMargin(bottom);
730   gPad->SetTopMargin(0.1);
731   gPad->SetBottomMargin(bottom);
732   gPad->SetRightMargin(0.2);
733   fhDcaXYvsDcaZnorm[1]->Draw("COLZ");
734
735   canvas3->SaveAs(Form("%s.eps", canvas3->GetName()));
736   canvas3->SaveAs(Form("%s.ps", canvas3->GetName()));
737
738   // -----
739
740   TCanvas* canvas4 = new TCanvas("Track_QA_Primary_4", "Track QA Primary 4", 1200, 500);
741   canvas4->Divide(3, 1);
742
743   canvas4->cd(1);
744   gPad->SetRightMargin(right);
745   gPad->SetLeftMargin(left);
746   gPad->SetTopMargin(top);
747   gPad->SetBottomMargin(bottom);
748   fhQA[kCutNSigmaToVertex][0]->SetStats(kFALSE);
749   fhQA[kCutNSigmaToVertex][0]->Draw();
750   fhQA[kCutNSigmaToVertex][1]->Draw("same");
751
752   canvas4->cd(2);
753   gPad->SetRightMargin(right);
754   gPad->SetLeftMargin(left);
755   gPad->SetTopMargin(top);
756   gPad->SetBottomMargin(bottom);
757   fhQA[kCutRequireSigmaToVertex][0]->SetStats(kFALSE);
758   fhQA[kCutRequireSigmaToVertex][0]->Draw();
759   fhQA[kCutRequireSigmaToVertex][1]->Draw("same");
760
761   canvas4->cd(3);
762   gPad->SetRightMargin(right);
763   gPad->SetLeftMargin(left);
764   gPad->SetTopMargin(top);
765   gPad->SetBottomMargin(bottom);
766   fhQA[kCutAcceptKinkDaughters][0]->SetStats(kFALSE);
767   fhQA[kCutAcceptKinkDaughters][0]->Draw();
768   fhQA[kCutAcceptKinkDaughters][1]->Draw("same");
769
770   canvas4->SaveAs(Form("%s.eps", canvas4->GetName()));
771   canvas4->SaveAs(Form("%s.ps", canvas4->GetName()));
772 }
773 //__________________________________________________________________________________
774 void AliCFTrackIsPrimaryCuts::AddQAHistograms(TList *qaList) {
775   //
776   // saves the histograms in a TList
777   //
778   DefineHistograms();
779
780   qaList->Add(fhCutStatistics);
781   qaList->Add(fhCutCorrelation);
782
783   for (Int_t j=0; j<kNStepQA; j++) {
784     qaList->Add(fhDcaXYvsDcaZ[j]);
785     qaList->Add(fhDcaXYvsDcaZnorm[j]);
786     for(Int_t i=0; i<kNHist; i++)
787         qaList->Add(fhQA[i][j]);
788   }
789 }