]> git.uio.no Git - u/mrichter/AliRoot.git/blob - CORRFW/AliCFTrackIsPrimaryCuts.cxx
Fix axis label
[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 // - min. and max. distance to main vertex in transverse plane (xy)
26 // - min. and max. longitudinal distance to main vertex (z)
27 // - min. and max. distance to main vertex as ellpise in xy - z plane
28 // - all above cuts on absolute values or in units of sigma (resolution)
29 // - min. and max. distance to main vertex in units of sigma (resolution)
30 // - max. transverse (xy) and longitudinal (z) impact parameter resolution
31 // - require that the dca calculation doesn't fail
32 // - accept or not accept daughter tracks of kink decays
33 //
34 // By default, the distance to 'vertex calculated from tracks' is used.
35 // Optionally the SPD (tracklet based) or TPC (TPC only tracks based) vertex
36 // can be used.
37 // Note: the distance to the TPC-vertex is already stored in the ESD,
38 // the distance to the SPD-vertex has to be re-calculated by propagating each
39 // track while executing this cut.
40 //
41 // The cut values for these cuts are set with the corresponding set functions.
42 // All cut classes provided by the correction framework are supposed to be
43 // added in the Analysis Framwork's class AliAnalysisFilter and applied by
44 // the filter via a loop.
45 //
46 // author: I. Kraus (Ingrid.Kraus@cern.ch)
47 // idea taken form
48 // AliESDtrackCuts writte by Jan Fiete Grosse-Oetringhaus and
49 // AliRsnDaughterCut class written by A. Pulvirenti.
50
51 #include <TCanvas.h>
52 #include <TDirectory.h>
53 #include <TH2.h>
54 #include <TBits.h>
55
56 #include <AliESDtrack.h>
57 #include <AliAODTrack.h>
58 #include <AliESDEvent.h>
59 #include <AliAODEvent.h>
60 #include <AliLog.h>
61 #include "AliCFTrackIsPrimaryCuts.h"
62
63 ClassImp(AliCFTrackIsPrimaryCuts)
64
65 //__________________________________________________________________________________
66 AliCFTrackIsPrimaryCuts::AliCFTrackIsPrimaryCuts() :
67   AliCFCutBase(),
68   fEvt(0x0),
69   fUseSPDvertex(0),
70   fUseTPCvertex(0),
71   fMinDCAToVertexXY(0),
72   fMinDCAToVertexZ(0),
73   fMaxDCAToVertexXY(0),
74   fMaxDCAToVertexZ(0),
75   fDCAToVertex2D(0),
76   fAbsDCAToVertex(0),
77   fNSigmaToVertexMin(0),
78   fNSigmaToVertexMax(0),
79   fSigmaDCAxy(0),
80   fSigmaDCAz(0),
81   fRequireSigmaToVertex(0),
82   fAODType(AliAODTrack::kUndef),
83   fAcceptKinkDaughters(0),
84   fhCutStatistics(0),
85   fhCutCorrelation(0),
86   fBitmap(0x0),
87   fhNBinsNSigma(0),
88   fhNBinsRequireSigma(0),
89   fhNBinsAcceptKink(0),
90   fhNBinsDcaXY(0),
91   fhNBinsDcaZ(0),
92   fhNBinsDcaXYnorm(0),
93   fhNBinsDcaZnorm(0),
94   fhNBinsSigmaDcaXY(0),
95   fhNBinsSigmaDcaZ(0),
96   fhBinLimNSigma(0x0),
97   fhBinLimRequireSigma(0x0),
98   fhBinLimAcceptKink(0x0),
99   fhBinLimDcaXY(0x0),
100   fhBinLimDcaZ(0x0),
101   fhBinLimDcaXYnorm(0x0),
102   fhBinLimDcaZnorm(0x0),
103   fhBinLimSigmaDcaXY(0x0),
104   fhBinLimSigmaDcaZ(0x0)
105 {
106   //
107   // Default constructor
108   //
109   Initialise();
110 }
111 //__________________________________________________________________________________
112 AliCFTrackIsPrimaryCuts::AliCFTrackIsPrimaryCuts(Char_t* name, Char_t* title) :
113   AliCFCutBase(name,title),
114   fEvt(0x0),
115   fUseSPDvertex(0),
116   fUseTPCvertex(0),
117   fMinDCAToVertexXY(0),
118   fMinDCAToVertexZ(0),
119   fMaxDCAToVertexXY(0),
120   fMaxDCAToVertexZ(0),
121   fDCAToVertex2D(0),
122   fAbsDCAToVertex(0),
123   fNSigmaToVertexMin(0),
124   fNSigmaToVertexMax(0),
125   fSigmaDCAxy(0),
126   fSigmaDCAz(0),
127   fRequireSigmaToVertex(0),
128   fAODType(AliAODTrack::kUndef),
129   fAcceptKinkDaughters(0),
130   fhCutStatistics(0),
131   fhCutCorrelation(0),
132   fBitmap(0x0),
133   fhNBinsNSigma(0),
134   fhNBinsRequireSigma(0),
135   fhNBinsAcceptKink(0),
136   fhNBinsDcaXY(0),
137   fhNBinsDcaZ(0),
138   fhNBinsDcaXYnorm(0),
139   fhNBinsDcaZnorm(0),
140   fhNBinsSigmaDcaXY(0),
141   fhNBinsSigmaDcaZ(0),
142   fhBinLimNSigma(0x0),
143   fhBinLimRequireSigma(0x0),
144   fhBinLimAcceptKink(0x0),
145   fhBinLimDcaXY(0x0),
146   fhBinLimDcaZ(0x0),
147   fhBinLimDcaXYnorm(0x0),
148   fhBinLimDcaZnorm(0x0),
149   fhBinLimSigmaDcaXY(0x0),
150   fhBinLimSigmaDcaZ(0x0)
151 {
152   //
153   // Constructor
154   //
155   Initialise();
156 }
157 //__________________________________________________________________________________
158 AliCFTrackIsPrimaryCuts::AliCFTrackIsPrimaryCuts(const AliCFTrackIsPrimaryCuts& c) :
159   AliCFCutBase(c),
160   fEvt(c.fEvt),
161   fUseSPDvertex(c.fUseSPDvertex),
162   fUseTPCvertex(c.fUseTPCvertex),
163   fMinDCAToVertexXY(c.fMinDCAToVertexXY),
164   fMinDCAToVertexZ(c.fMinDCAToVertexZ),
165   fMaxDCAToVertexXY(c.fMaxDCAToVertexXY),
166   fMaxDCAToVertexZ(c.fMaxDCAToVertexZ),
167   fDCAToVertex2D(c.fDCAToVertex2D),
168   fAbsDCAToVertex(c.fAbsDCAToVertex),
169   fNSigmaToVertexMin(c.fNSigmaToVertexMin),
170   fNSigmaToVertexMax(c.fNSigmaToVertexMax),
171   fSigmaDCAxy(c.fSigmaDCAxy),
172   fSigmaDCAz(c.fSigmaDCAz),
173   fRequireSigmaToVertex(c.fRequireSigmaToVertex),
174   fAODType(c.fAODType),
175   fAcceptKinkDaughters(c.fAcceptKinkDaughters),
176   fhCutStatistics(c.fhCutStatistics),
177   fhCutCorrelation(c.fhCutCorrelation),
178   fBitmap(c.fBitmap),
179   fhNBinsNSigma(c.fhNBinsNSigma),
180   fhNBinsRequireSigma(c.fhNBinsRequireSigma),
181   fhNBinsAcceptKink(c.fhNBinsAcceptKink),
182   fhNBinsDcaXY(c.fhNBinsDcaXY),
183   fhNBinsDcaZ(c.fhNBinsDcaZ),
184   fhNBinsDcaXYnorm(c.fhNBinsDcaXYnorm),
185   fhNBinsDcaZnorm(c.fhNBinsDcaZnorm),
186   fhNBinsSigmaDcaXY(c.fhNBinsSigmaDcaXY),
187   fhNBinsSigmaDcaZ(c.fhNBinsSigmaDcaZ),
188   fhBinLimNSigma(c.fhBinLimNSigma),
189   fhBinLimRequireSigma(c.fhBinLimRequireSigma),
190   fhBinLimAcceptKink(c.fhBinLimAcceptKink),
191   fhBinLimDcaXY(c.fhBinLimDcaXY),
192   fhBinLimDcaZ(c.fhBinLimDcaZ),
193   fhBinLimDcaXYnorm(c.fhBinLimDcaXYnorm),
194   fhBinLimDcaZnorm(c.fhBinLimDcaZnorm),
195   fhBinLimSigmaDcaXY(c.fhBinLimSigmaDcaXY),
196   fhBinLimSigmaDcaZ(c.fhBinLimSigmaDcaZ)
197 {
198   //
199   // copy constructor
200   //
201   ((AliCFTrackIsPrimaryCuts &) c).Copy(*this);
202 }
203 //__________________________________________________________________________________
204 AliCFTrackIsPrimaryCuts& AliCFTrackIsPrimaryCuts::operator=(const AliCFTrackIsPrimaryCuts& c)
205 {
206   //
207   // Assignment operator
208   //
209   if (this != &c) {
210     AliCFCutBase::operator=(c) ;
211     fEvt = c.fEvt;
212     fUseSPDvertex = c.fUseSPDvertex;
213     fUseTPCvertex = c.fUseTPCvertex;
214     fMinDCAToVertexXY = c.fMinDCAToVertexXY;
215     fMinDCAToVertexZ = c.fMinDCAToVertexZ;
216     fMaxDCAToVertexXY = c.fMaxDCAToVertexXY;
217     fMaxDCAToVertexZ = c.fMaxDCAToVertexZ;
218     fDCAToVertex2D = c.fDCAToVertex2D;
219     fAbsDCAToVertex = c.fAbsDCAToVertex;
220     fNSigmaToVertexMin = c.fNSigmaToVertexMin ;
221     fNSigmaToVertexMax = c.fNSigmaToVertexMax ;
222     fSigmaDCAxy = c.fSigmaDCAxy ;
223     fSigmaDCAz = c.fSigmaDCAz ;
224     fRequireSigmaToVertex = c.fRequireSigmaToVertex ;
225     fAODType = c.fAODType ;
226     fAcceptKinkDaughters = c.fAcceptKinkDaughters ;
227     fhCutStatistics = c.fhCutStatistics ;
228     fhCutCorrelation = c.fhCutCorrelation ;
229     fBitmap =  c.fBitmap;
230     fhNBinsNSigma = c.fhNBinsNSigma;
231     fhNBinsRequireSigma = c.fhNBinsRequireSigma;
232     fhNBinsAcceptKink = c.fhNBinsAcceptKink;
233     fhNBinsDcaXY = c.fhNBinsDcaXY;
234     fhNBinsDcaZ = c.fhNBinsDcaZ;
235     fhNBinsDcaXYnorm = c.fhNBinsDcaXYnorm;
236     fhNBinsDcaZnorm = c.fhNBinsDcaZnorm;
237     fhNBinsSigmaDcaXY = c.fhNBinsSigmaDcaXY;
238     fhNBinsSigmaDcaZ = c.fhNBinsSigmaDcaZ;
239     fhBinLimNSigma = c.fhBinLimNSigma;
240     fhBinLimRequireSigma = c.fhBinLimRequireSigma;
241     fhBinLimAcceptKink = c.fhBinLimAcceptKink;
242     fhBinLimDcaXY = c.fhBinLimDcaXY;
243     fhBinLimDcaZ = c.fhBinLimDcaZ;
244     fhBinLimDcaXYnorm = c.fhBinLimDcaXYnorm;
245     fhBinLimDcaZnorm = c.fhBinLimDcaZnorm;
246     fhBinLimSigmaDcaXY = c.fhBinLimSigmaDcaXY;
247     fhBinLimSigmaDcaZ = c.fhBinLimSigmaDcaZ;
248
249     for (Int_t j=0; j<6; j++) fDCA[j] = c.fDCA[j];
250     for (Int_t j=0; j<c.kNStepQA; j++){
251       if(c.fhDcaXYvsDcaZ[j]) fhDcaXYvsDcaZ[j] = (TH2F*)c.fhDcaXYvsDcaZ[j]->Clone();
252       for (Int_t i=0; i<c.kNHist; i++){
253         if(c.fhQA[i][j]) fhQA[i][j] = (TH1F*)c.fhQA[i][j]->Clone();
254       }
255     }
256     ((AliCFTrackIsPrimaryCuts &) c).Copy(*this);
257  }
258   return *this;
259 }
260 //__________________________________________________________________________________
261 AliCFTrackIsPrimaryCuts::~AliCFTrackIsPrimaryCuts()
262 {
263   //
264   // destructor
265   //
266   if (fhCutStatistics)                  delete fhCutStatistics;
267   if (fhCutCorrelation)                 delete fhCutCorrelation;
268
269   for (Int_t j=0; j<kNStepQA; j++){
270     if(fhDcaXYvsDcaZ[j])        delete fhDcaXYvsDcaZ[j];
271     for (Int_t i=0; i<kNHist; i++)
272       if(fhQA[i][j])            delete fhQA[i][j];
273   }
274   if(fEvt)                      delete fEvt;
275   if(fBitmap)                   delete fBitmap;
276   if(fhBinLimNSigma)            delete fhBinLimNSigma;
277   if(fhBinLimRequireSigma)      delete fhBinLimRequireSigma;
278   if(fhBinLimAcceptKink)        delete fhBinLimAcceptKink;
279   if(fhBinLimDcaXY)             delete fhBinLimDcaXY;
280   if(fhBinLimDcaZ)              delete fhBinLimDcaZ;
281   if(fhBinLimDcaXYnorm)         delete fhBinLimDcaXYnorm;
282   if(fhBinLimDcaZnorm)          delete fhBinLimDcaZnorm;
283   if(fhBinLimSigmaDcaXY)        delete fhBinLimSigmaDcaXY;
284   if(fhBinLimSigmaDcaZ)         delete fhBinLimSigmaDcaZ;
285 }
286 //__________________________________________________________________________________
287 void AliCFTrackIsPrimaryCuts::Initialise()
288 {
289   //
290   // sets everything to zero
291   //
292   fUseSPDvertex = 0;
293   fUseTPCvertex = 0;
294   fMinDCAToVertexXY = 0;
295   fMinDCAToVertexZ = 0;
296   fMaxDCAToVertexXY = 0;
297   fMaxDCAToVertexZ = 0;
298   fDCAToVertex2D = 0;
299   fAbsDCAToVertex = 0;
300   fNSigmaToVertexMin = 0;
301   fNSigmaToVertexMax = 0;
302   fSigmaDCAxy = 0;
303   fSigmaDCAz = 0;
304   fRequireSigmaToVertex = 0;
305   fAcceptKinkDaughters = 0;
306   fAODType = AliAODTrack::kUndef;
307
308   SetMinDCAToVertexXY();
309   SetMinDCAToVertexZ();
310   SetMaxDCAToVertexXY();
311   SetMaxDCAToVertexZ();
312   SetDCAToVertex2D();
313   SetAbsDCAToVertex();
314   SetMinNSigmaToVertex();
315   SetMaxNSigmaToVertex();
316   SetMaxSigmaDCAxy();
317   SetMaxSigmaDCAz();
318   SetRequireSigmaToVertex();
319   SetAcceptKinkDaughters();
320   SetAODType();
321
322   for (Int_t j=0; j<6; j++) fDCA[j] = 0.;
323   for (Int_t j=0; j<kNStepQA; j++)  {
324     fhDcaXYvsDcaZ[j] = 0x0;
325     for (Int_t i=0; i<kNHist; i++)
326       fhQA[i][j] = 0x0;
327   }
328   fhCutStatistics = 0;
329   fhCutCorrelation = 0;
330   fBitmap=new TBits(0);
331
332   //set default bining for QA histograms
333   SetHistogramBins(kCutNSigmaToVertex,100,0.,10.);
334   SetHistogramBins(kCutRequireSigmaToVertex,5,-0.75,1.75);
335   SetHistogramBins(kCutAcceptKinkDaughters,5,-0.75,1.75);
336   SetHistogramBins(kDcaXY,500,-10.,10.);
337   SetHistogramBins(kDcaZ,500,-10.,10.);
338   SetHistogramBins(kDcaXYnorm,500,-10.,10.);
339   SetHistogramBins(kDcaZnorm,500,-10.,10.);
340   SetHistogramBins(kSigmaDcaXY,500,-0.1,0.9);
341   SetHistogramBins(kSigmaDcaZ,500,-0.1,0.9);
342 }
343 //__________________________________________________________________________________
344 void AliCFTrackIsPrimaryCuts::Copy(TObject &c) const
345 {
346   //
347   // Copy function
348   //
349   AliCFTrackIsPrimaryCuts& target = (AliCFTrackIsPrimaryCuts &) c;
350
351   target.Initialise();
352
353   if (fhCutStatistics)  target.fhCutStatistics = (TH1F*) fhCutStatistics->Clone();
354   if (fhCutCorrelation) target.fhCutCorrelation = (TH2F*) fhCutCorrelation->Clone();
355
356   for (Int_t j=0; j<6; j++) target.fDCA[j] = fDCA[j];
357   for (Int_t j=0; j<kNStepQA; j++){
358     if(fhDcaXYvsDcaZ[j]) target.fhDcaXYvsDcaZ[j] = (TH2F*)fhDcaXYvsDcaZ[j]->Clone();
359     for (Int_t i=0; i<kNHist; i++)
360       if(fhQA[i][j]) target.fhQA[i][j] = (TH1F*)fhQA[i][j]->Clone();
361   }
362   TNamed::Copy(c);
363 }
364 //__________________________________________________________________________________
365 void AliCFTrackIsPrimaryCuts::SetRecEventInfo(const TObject* evt) {
366   //
367   // Sets pointer to event information (AliESDEvent or AliAODEvent)
368   //
369   if (!evt) {
370     AliError("Pointer to AliVEvent !");
371     return;
372   }
373   TString className(evt->ClassName());
374   if (! (className.CompareTo("AliESDEvent")==0 || className.CompareTo("AliAODEvent")==0)) {
375     AliError("argument must point to an AliESDEvent or AliAODEvent !");
376     return ;
377   }
378   fEvt = (AliVEvent*) evt;
379 }
380 //__________________________________________________________________________________
381 void AliCFTrackIsPrimaryCuts::UseSPDvertex(Bool_t b) {
382   fUseSPDvertex = b;
383   if(fUseTPCvertex && fUseSPDvertex) {
384         fUseSPDvertex = kFALSE;
385         AliError("SPD and TPC vertex chosen. TPC vertex is preferred.");
386   }
387 }
388 //__________________________________________________________________________________
389 void AliCFTrackIsPrimaryCuts::UseTPCvertex(Bool_t b) {
390   fUseTPCvertex = b;
391   if(fUseTPCvertex) fUseSPDvertex = kFALSE;
392 }
393 //__________________________________________________________________________________
394 void AliCFTrackIsPrimaryCuts::GetDCA(AliESDtrack* esdTrack)
395 {
396   if (!esdTrack) return;
397
398   Float_t b[2] = {0.,0.};
399   Float_t bCov[3] = {0.,0.,0.};
400   if(!fUseSPDvertex && !fUseTPCvertex)  esdTrack->GetImpactParameters(b,bCov);
401   if( fUseTPCvertex)                    esdTrack->GetImpactParametersTPC(b,bCov);
402
403   if( fUseSPDvertex) {
404         if (!fEvt) return;
405         AliESDEvent * evt = 0x0 ; 
406         evt = dynamic_cast<AliESDEvent*>(fEvt);
407         if (!evt) {
408           AliError("event not found");
409           return;
410         }
411         const AliESDVertex *vtx = evt->GetVertex();
412         const Double_t Bz = evt->GetMagneticField();
413         AliExternalTrackParam *cParam = 0;
414         Bool_t success = esdTrack->RelateToVertex(vtx, Bz, kVeryBig, cParam);
415         if (success) esdTrack->GetImpactParameters(b,bCov);
416   }
417
418   if (bCov[0]<=0 || bCov[2]<=0) {
419       bCov[0]=0; bCov[2]=0;
420   }
421   fDCA[0] = b[0]; // impact parameter xy
422   fDCA[1] = b[1]; // impact parameter z
423   fDCA[2] = TMath::Sqrt(bCov[0]); // resolution xy
424   fDCA[3] = TMath::Sqrt(bCov[2]); // resolution z
425
426   if (!fAbsDCAToVertex) {
427         if (fDCA[2] > 0) fDCA[0] = fDCA[0]/fDCA[2]; // normalised impact parameter xy
428         if (fDCA[3] > 0) fDCA[1] = fDCA[1]/fDCA[3]; // normalised impact parameter z
429   }
430
431   // get n_sigma
432   if(!fUseSPDvertex && !fUseTPCvertex)
433         fDCA[5] = AliESDtrackCuts::GetSigmaToVertex(esdTrack);
434
435   if(fUseTPCvertex) {
436         fDCA[5] = -1;
437         if (fDCA[2]==0 || fDCA[3]==0)
438           return;
439         fDCA[5] = 1000.;
440         Float_t d = TMath::Sqrt(TMath::Power(b[0]/fDCA[2],2) + TMath::Power(b[1]/fDCA[3],2));
441         if (TMath::Exp(-d * d / 2) < 1e-15)
442           return;
443         fDCA[5] = TMath::ErfInverse(1 - TMath::Exp(-d * d / 2)) * TMath::Sqrt(2);
444   }
445   return;
446 }
447 //__________________________________________________________________________________
448 void AliCFTrackIsPrimaryCuts::GetDCA(AliAODTrack* aodTrack)
449 {
450   if (!aodTrack) return;
451
452   Double_t p[3] = {0.,0.,0.};
453   Double_t cov[21] = {0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
454
455   aodTrack->XYZAtDCA(p); // position at DCA
456   aodTrack->GetCovarianceXYZPxPyPz(cov);
457
458
459   fDCA[5] = -1; // n_sigma
460   if (p[0]==-999. || p[1]==-999. ||  p[2]==-999.) {
461     AliError("dca info not available !");
462     fDCA[0] = -999.; // impact parameter xy
463     fDCA[1] = -999.; // impact parameter z
464     return;
465   }
466
467   AliAODEvent * evt = 0x0;
468   evt = dynamic_cast<AliAODEvent*>(fEvt);
469   if (!evt) return;
470
471   // primary vertex is the "best": tracks, SPD or TPC vertex
472   AliAODVertex * primaryVertex = evt->GetVertex(0);
473   // dca = track postion - primary vertex position
474   p[0] = p[0] - primaryVertex->GetX();
475   p[1] = p[1] - primaryVertex->GetY();
476   p[2] = p[2] - primaryVertex->GetZ();
477   // calculate dca in transverse plane
478   Float_t b[2] = {0.,0.};
479   b[0] = TMath::Sqrt(p[0]*p[0]+p[1]*p[1]);
480   b[1] = p[2];
481   // resolution
482   Double_t bCov[3] = {0.,0.,0.};
483   // how to calculate the resoultion in the transverse plane ?
484   bCov[0] = 0.; // to do: calculate cov in transverse plane
485   bCov[2] = 0.; // from cov in x and y, need to know correlation
486
487   if (bCov[0]<=0 || bCov[2]<=0) {
488       bCov[0]=0; bCov[2]=0;
489   }
490   fDCA[0] = b[0]; // impact parameter xy
491   fDCA[1] = b[1]; // impact parameter z
492   fDCA[2] = TMath::Sqrt(bCov[0]); // resolution xy
493   fDCA[3] = TMath::Sqrt(bCov[2]); // resolution z
494
495   if (!fAbsDCAToVertex) {
496         AliError("resolution of impact parameter not available, use absolute dca cut instead !");
497         if (fDCA[2] > 0) fDCA[0] = fDCA[0]/fDCA[2]; // normalised impact parameter xy
498         if (fDCA[3] > 0) fDCA[1] = fDCA[1]/fDCA[3]; // normalised impact parameter z
499   }
500
501   // get n_sigma
502   if (fDCA[2]==0 || fDCA[3]==0)
503         return;
504   fDCA[5] = 1000.;
505   Float_t d = TMath::Sqrt(TMath::Power(b[0]/fDCA[2],2) + TMath::Power(b[1]/fDCA[3],2));
506   if (TMath::Exp(-d * d / 2) < 1e-15)
507         return;
508   fDCA[5] = TMath::ErfInverse(1 - TMath::Exp(-d * d / 2)) * TMath::Sqrt(2);
509
510   return;
511 }
512 //__________________________________________________________________________________
513 void AliCFTrackIsPrimaryCuts::SelectionBitMap(TObject* obj)
514 {
515   //
516   // test if the track passes the single cuts
517   // and store the information in a bitmap
518   //
519
520   // bitmap stores the decision of each single cut
521   for(Int_t i=0; i<kNCuts; i++)fBitmap->SetBitNumber(i,kFALSE);
522
523   // check TObject and cast into ESDtrack
524   if (!obj) return;
525   if (!obj->InheritsFrom("AliVParticle")) {
526     AliError("object must derived from AliVParticle !");
527     return;
528   }
529   
530   AliESDtrack * esdTrack = dynamic_cast<AliESDtrack*>(obj);
531   AliAODTrack * aodTrack = dynamic_cast<AliAODTrack*>(obj);
532
533   if (!(esdTrack || aodTrack)) {
534     AliError("object must be an ESDtrack or an AODtrack !");
535     return;
536   }
537
538   Bool_t isESDTrack = kFALSE;
539   Bool_t isAODTrack = kFALSE;
540
541   if (esdTrack) isESDTrack = strcmp(obj->ClassName(),"AliESDtrack") == 0 ? kTRUE : kFALSE ;
542   if (aodTrack) isAODTrack = strcmp(obj->ClassName(),"AliAODTrack") == 0 ? kTRUE : kFALSE ;
543
544   // get the track to vertex parameter for ESD track
545   if (isESDTrack) GetDCA(esdTrack);
546   if (isAODTrack) GetDCA(aodTrack);
547
548   // check whether dca info is filled
549   Bool_t dcaInfo = 0;
550   if (fDCA[0]>-990. && fDCA[1]>-990.) dcaInfo = 1;
551
552   Float_t bxy = 0, bz = 0;
553   bxy = TMath::Abs(fDCA[0]);
554   bz  = TMath::Abs(fDCA[1]);
555
556   Float_t b2Dmin = 0, b2Dmax = 0;
557   if (fMinDCAToVertexXY>0 && fMinDCAToVertexZ>0)
558     b2Dmin = fDCA[0]*fDCA[0]/fMinDCAToVertexXY/fMinDCAToVertexXY + fDCA[1]*fDCA[1]/fMinDCAToVertexZ/fMinDCAToVertexZ;
559   if (fMaxDCAToVertexXY>0 && fMaxDCAToVertexZ>0) 
560     b2Dmax = fDCA[0]*fDCA[0]/fMaxDCAToVertexXY/fMaxDCAToVertexXY + fDCA[1]*fDCA[1]/fMaxDCAToVertexZ/fMaxDCAToVertexZ;
561
562
563   // fill the bitmap
564   Int_t iCutBit = 0;
565
566   if (!dcaInfo || fDCAToVertex2D || (!fDCAToVertex2D && bxy >= fMinDCAToVertexXY && bxy <= fMaxDCAToVertexXY))
567         fBitmap->SetBitNumber(iCutBit,kTRUE);
568   iCutBit++;
569
570   if (!dcaInfo || fDCAToVertex2D || (!fDCAToVertex2D && bz  >= fMinDCAToVertexZ && bz  <= fMaxDCAToVertexZ))
571         fBitmap->SetBitNumber(iCutBit,kTRUE);
572   iCutBit++;
573
574   if (!dcaInfo || !fDCAToVertex2D || (fDCAToVertex2D && TMath::Sqrt(b2Dmin) > 1  && TMath::Sqrt(b2Dmax) < 1))
575       fBitmap->SetBitNumber(iCutBit,kTRUE);
576   iCutBit++;
577
578   if (!dcaInfo || (fDCA[5] >= fNSigmaToVertexMin && fDCA[5] <= fNSigmaToVertexMax))
579       fBitmap->SetBitNumber(iCutBit,kTRUE);
580   iCutBit++;
581
582   if (!dcaInfo || fDCA[2] < fSigmaDCAxy)
583       fBitmap->SetBitNumber(iCutBit,kTRUE);
584   iCutBit++;
585
586   if (!dcaInfo || fDCA[3] < fSigmaDCAz)
587       fBitmap->SetBitNumber(iCutBit,kTRUE);
588   iCutBit++;
589
590   if (!dcaInfo || !fRequireSigmaToVertex || (fDCA[5]>=0 && fRequireSigmaToVertex))
591       fBitmap->SetBitNumber(iCutBit,kTRUE);
592   iCutBit++;
593
594   if (!dcaInfo || fAcceptKinkDaughters || (!fAcceptKinkDaughters && esdTrack->GetKinkIndex(0)<=0))
595       fBitmap->SetBitNumber(iCutBit,kTRUE);
596   iCutBit++;
597
598   if (isAODTrack) {
599     if (fAODType==AliAODTrack::kUndef || fAODType == aodTrack->GetType()) {
600       fBitmap->SetBitNumber(iCutBit,kTRUE);
601     }
602   }
603   else fBitmap->SetBitNumber(iCutBit,kTRUE);
604
605   return;
606 }
607 //__________________________________________________________________________________
608 Bool_t AliCFTrackIsPrimaryCuts::IsSelected(TObject* obj) {
609   //
610   // loops over decisions of single cuts and returns if the track is accepted
611   //
612   SelectionBitMap(obj);
613
614   if (fIsQAOn) FillHistograms(obj,0);
615   Bool_t isSelected = kTRUE;
616
617   for (UInt_t icut=0; icut<fBitmap->GetNbits();icut++) {
618     if(!fBitmap->TestBitNumber(icut)) {
619         isSelected = kFALSE;
620         break;
621     }
622   }
623
624   if (!isSelected) return kFALSE ;
625   if (fIsQAOn) FillHistograms(obj,1);
626   return kTRUE;
627 }
628 //__________________________________________________________________________________
629 void AliCFTrackIsPrimaryCuts::SetHistogramBins(Int_t index, Int_t nbins, Double_t *bins)
630 {
631   //
632   // variable bin size
633   //
634   if(!fIsQAOn) return;
635
636   switch(index){
637   case kCutNSigmaToVertex:
638     fhNBinsNSigma=nbins+1;
639     fhBinLimNSigma=new Double_t[nbins+1];
640     for(Int_t i=0;i<nbins+1;i++)fhBinLimNSigma[i]=bins[i];
641     break;
642
643   case kCutRequireSigmaToVertex:
644     fhNBinsRequireSigma=nbins+1;
645     fhBinLimRequireSigma=new Double_t[nbins+1];
646     for(Int_t i=0;i<nbins+1;i++)fhBinLimRequireSigma[i]=bins[i];
647     break;
648
649   case kCutAcceptKinkDaughters:
650     fhNBinsAcceptKink=nbins+1;
651     fhBinLimAcceptKink=new Double_t[nbins+1];
652     for(Int_t i=0;i<nbins+1;i++)fhBinLimAcceptKink[i]=bins[i];
653     break;
654
655   case kDcaXY:
656     fhNBinsDcaXY=nbins+1;
657     fhBinLimDcaXY=new Double_t[nbins+1];
658     for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaXY[i]=bins[i];
659     break;
660
661   case kDcaZ:
662     fhNBinsDcaZ=nbins+1;
663     fhBinLimDcaZ=new Double_t[nbins+1];
664     for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaZ[i]=bins[i];
665     break;
666
667   case kDcaXYnorm:
668     fhNBinsDcaXYnorm=nbins+1;
669     fhBinLimDcaXYnorm=new Double_t[nbins+1];
670     for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaXYnorm[i]=bins[i];
671     break;
672
673   case kDcaZnorm:
674     fhNBinsDcaZnorm=nbins+1;
675     fhBinLimDcaZnorm=new Double_t[nbins+1];
676     for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaZnorm[i]=bins[i];
677     break;
678
679   case kSigmaDcaXY:
680     fhNBinsSigmaDcaXY=nbins+1;
681     fhBinLimSigmaDcaXY=new Double_t[nbins+1];
682     for(Int_t i=0;i<nbins+1;i++)fhBinLimSigmaDcaXY[i]=bins[i];
683     break;
684
685   case kSigmaDcaZ:
686     fhNBinsSigmaDcaZ=nbins+1;
687     fhBinLimSigmaDcaZ=new Double_t[nbins+1];
688     for(Int_t i=0;i<nbins+1;i++)fhBinLimSigmaDcaZ[i]=bins[i];
689     break;
690   }
691 }
692 //__________________________________________________________________________________
693 void AliCFTrackIsPrimaryCuts::SetHistogramBins(Int_t index, Int_t nbins, Double_t xmin, Double_t xmax)
694 {
695   //
696   // fixed bin size
697   //
698   switch(index){
699   case kCutNSigmaToVertex:
700     fhNBinsNSigma=nbins+1;
701     fhBinLimNSigma=new Double_t[nbins+1];
702     for(Int_t i=0;i<nbins+1;i++)fhBinLimNSigma[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
703     break;
704
705   case kCutRequireSigmaToVertex:
706     fhNBinsRequireSigma=nbins+1;
707     fhBinLimRequireSigma=new Double_t[nbins+1];
708     for(Int_t i=0;i<nbins+1;i++)fhBinLimRequireSigma[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
709     break;
710
711   case kCutAcceptKinkDaughters:
712     fhNBinsAcceptKink=nbins+1;
713     fhBinLimAcceptKink=new Double_t[nbins+1];
714     for(Int_t i=0;i<nbins+1;i++)fhBinLimAcceptKink[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
715     break;
716
717   case kDcaXY:
718     fhNBinsDcaXY=nbins+1;
719     fhBinLimDcaXY=new Double_t[nbins+1];
720     for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaXY[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
721     break;
722
723   case kDcaZ:
724     fhNBinsDcaZ=nbins+1;
725     fhBinLimDcaZ=new Double_t[nbins+1];
726     for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaZ[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
727     break;
728
729   case kDcaXYnorm:
730     fhNBinsDcaXYnorm=nbins+1;
731     fhBinLimDcaXYnorm=new Double_t[nbins+1];
732     for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaXYnorm[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
733     break;
734
735   case kDcaZnorm:
736     fhNBinsDcaZnorm=nbins+1;
737     fhBinLimDcaZnorm=new Double_t[nbins+1];
738     for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaZnorm[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
739     break;
740
741   case kSigmaDcaXY:
742     fhNBinsSigmaDcaXY=nbins+1;
743     fhBinLimSigmaDcaXY=new Double_t[nbins+1];
744     for(Int_t i=0;i<nbins+1;i++)fhBinLimSigmaDcaXY[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
745     break;
746
747   case kSigmaDcaZ:
748     fhNBinsSigmaDcaZ=nbins+1;
749     fhBinLimSigmaDcaZ=new Double_t[nbins+1];
750     for(Int_t i=0;i<nbins+1;i++)fhBinLimSigmaDcaZ[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
751     break;
752   }
753 }
754 //__________________________________________________________________________________
755  void AliCFTrackIsPrimaryCuts::DefineHistograms() {
756   //
757   // histograms for cut variables, cut statistics and cut correlations
758   //
759   Int_t color = 2;
760
761   // book cut statistics and cut correlation histograms
762   fhCutStatistics = new TH1F(Form("%s_cut_statistics",GetName()), Form("%s cut statistics",GetName()), kNCuts,0.5,kNCuts+0.5);
763   fhCutStatistics->SetLineWidth(2);
764   fhCutStatistics->GetXaxis()->SetBinLabel(1,"dca xy");
765   fhCutStatistics->GetXaxis()->SetBinLabel(2,"dca z");
766   fhCutStatistics->GetXaxis()->SetBinLabel(3,"dca ellipse");
767   fhCutStatistics->GetXaxis()->SetBinLabel(4,"n dca");
768   fhCutStatistics->GetXaxis()->SetBinLabel(5,"sigma dca xy");
769   fhCutStatistics->GetXaxis()->SetBinLabel(6,"sigma dca z");
770   fhCutStatistics->GetXaxis()->SetBinLabel(7,"require dca");
771   fhCutStatistics->GetXaxis()->SetBinLabel(8,"kink daughter");
772   fhCutStatistics->GetXaxis()->SetBinLabel(9,"AOD type");
773
774   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);
775   fhCutCorrelation->SetLineWidth(2);
776   fhCutCorrelation->GetXaxis()->SetBinLabel(1,"dca xy");
777   fhCutCorrelation->GetXaxis()->SetBinLabel(2,"dca z");
778   fhCutCorrelation->GetXaxis()->SetBinLabel(3,"dca ellipse");
779   fhCutCorrelation->GetXaxis()->SetBinLabel(4,"n dca");
780   fhCutCorrelation->GetXaxis()->SetBinLabel(5,"sigma dca xy");
781   fhCutCorrelation->GetXaxis()->SetBinLabel(6,"sigma dca z");
782   fhCutCorrelation->GetXaxis()->SetBinLabel(7,"require dca");
783   fhCutCorrelation->GetXaxis()->SetBinLabel(8,"kink daughter");
784   fhCutCorrelation->GetXaxis()->SetBinLabel(9,"AOD type");
785
786   fhCutCorrelation->GetYaxis()->SetBinLabel(1,"dca xy");
787   fhCutCorrelation->GetYaxis()->SetBinLabel(2,"dca z");
788   fhCutCorrelation->GetYaxis()->SetBinLabel(3,"dca ellipse");
789   fhCutCorrelation->GetYaxis()->SetBinLabel(4,"n dca");
790   fhCutCorrelation->GetYaxis()->SetBinLabel(5,"sigma dca xy");
791   fhCutCorrelation->GetYaxis()->SetBinLabel(6,"sigma dca z");
792   fhCutCorrelation->GetYaxis()->SetBinLabel(7,"require dca");
793   fhCutCorrelation->GetYaxis()->SetBinLabel(8,"kink daughter");
794   fhCutCorrelation->GetYaxis()->SetBinLabel(9,"AOD type");
795
796   // book QA histograms
797   Char_t str[5];
798   for (Int_t i=0; i<kNStepQA; i++) {
799     if (i==0) snprintf(str,5," ");
800     else snprintf(str,5,"_cut");
801
802     fhDcaXYvsDcaZ[i]            = new  TH2F(Form("%s_dcaXYvsDcaZ%s",GetName(),str),"",200,-10,10,200,-10,10);
803     fhQA[kCutNSigmaToVertex][i] = new TH1F(Form("%s_nSigmaToVertex%s",GetName(),str),"",fhNBinsNSigma-1,fhBinLimNSigma);
804     fhQA[kCutRequireSigmaToVertex][i] = new TH1F(Form("%s_requireSigmaToVertex%s",GetName(),str),"",fhNBinsRequireSigma-1,fhBinLimRequireSigma);
805     fhQA[kCutAcceptKinkDaughters][i] = new TH1F(Form("%s_acceptKinkDaughters%s",GetName(),str),"",fhNBinsAcceptKink-1,fhBinLimAcceptKink);
806     fhQA[kDcaXY][i]             = new TH1F(Form("%s_dcaXY%s",GetName(),str),"",fhNBinsDcaXY-1,fhBinLimDcaXY);
807     fhQA[kDcaZ][i]              = new TH1F(Form("%s_dcaZ%s",GetName(),str),"",fhNBinsDcaZ-1,fhBinLimDcaZ);
808     fhQA[kDcaXYnorm][i]         = new TH1F(Form("%s_dcaXYnorm%s",GetName(),str),"",fhNBinsDcaXYnorm-1,fhBinLimDcaXYnorm);
809     fhQA[kDcaZnorm][i]          = new TH1F(Form("%s_dcaZnorm%s",GetName(),str),"",fhNBinsDcaZnorm-1,fhBinLimDcaZnorm);
810     fhQA[kSigmaDcaXY][i]        = new TH1F(Form("%s_sigmaDcaXY%s",GetName(),str),"",fhNBinsSigmaDcaXY-1,fhBinLimSigmaDcaXY);
811     fhQA[kSigmaDcaZ][i]         = new TH1F(Form("%s_sigmaDcaZ%s",GetName(),str),"",fhNBinsSigmaDcaZ-1,fhBinLimSigmaDcaZ);
812
813     fhDcaXYvsDcaZ[i]->SetXTitle("impact par. d_{z}");
814     fhDcaXYvsDcaZ[i]->SetYTitle("impact par. d_{xy}");
815
816     fhQA[kCutNSigmaToVertex][i]->SetXTitle("n #sigma to vertex");
817     fhQA[kCutRequireSigmaToVertex][i]->SetXTitle("require #sigma to vertex");
818     fhQA[kCutAcceptKinkDaughters][i]->SetXTitle("accept kink daughters");
819     fhQA[kDcaXY][i]->SetXTitle("impact par. d_{xy}");
820     fhQA[kDcaZ][i]->SetXTitle("impact par. d_{z}");
821     fhQA[kDcaXYnorm][i]->SetXTitle("norm. impact par. d_{xy} / #sigma_{xy}");
822     fhQA[kDcaZnorm][i]->SetXTitle("norm. impact par. d_{z} / #sigma_{z}");
823     fhQA[kSigmaDcaXY][i]->SetXTitle("impact par. resolution #sigma_{xy}");
824     fhQA[kSigmaDcaZ][i]->SetXTitle("impact par. resolution #sigma_{z}");
825   }
826   for(Int_t i=0; i<kNHist; i++) fhQA[i][1]->SetLineColor(color);
827 }
828 //__________________________________________________________________________________
829 void AliCFTrackIsPrimaryCuts::FillHistograms(TObject* obj, Bool_t f)
830 {
831   //
832   // fill the QA histograms
833   //
834
835   if (!obj) return;
836   Bool_t isESDTrack = strcmp(obj->ClassName(),"AliESDtrack") == 0 ? kTRUE : kFALSE ;
837   Bool_t isAODTrack = strcmp(obj->ClassName(),"AliAODTrack") == 0 ? kTRUE : kFALSE ;
838
839   AliESDtrack * esdTrack = 0x0 ;
840   AliAODTrack * aodTrack = 0x0 ; 
841   if (isESDTrack) esdTrack = dynamic_cast<AliESDtrack*>(obj);
842   if (isAODTrack) aodTrack = dynamic_cast<AliAODTrack*>(obj);
843
844   // f = 0: fill histograms before cuts
845   // f = 1: fill histograms after cuts
846
847   // get the track to vertex parameter for ESD track
848   if (esdTrack) {
849
850     fhQA[kDcaZ][f]->Fill(fDCA[1]);
851     fhQA[kDcaXY][f]->Fill(fDCA[0]);
852     fhDcaXYvsDcaZ[f]->Fill(fDCA[1],fDCA[0]);
853     fhQA[kSigmaDcaXY][f]->Fill(fDCA[2]);
854     fhQA[kSigmaDcaZ][f]->Fill(fDCA[3]);
855 // // // // // // //  delete histograms
856       fhQA[kDcaZnorm][f]->Fill(fDCA[1]);
857       fhQA[kDcaXYnorm][f]->Fill(fDCA[0]);
858
859     fhQA[kCutNSigmaToVertex][f]->Fill(fDCA[5]);
860     if (fDCA[5]<0 && fRequireSigmaToVertex) 
861       fhQA[kCutRequireSigmaToVertex][f]->Fill(0.);
862     else 
863       fhQA[kCutRequireSigmaToVertex][f]->Fill(1.);
864
865     if (!fAcceptKinkDaughters && esdTrack->GetKinkIndex(0)>0)
866       fhQA[kCutAcceptKinkDaughters][f]->Fill(0.);
867     else 
868       fhQA[kCutAcceptKinkDaughters][f]->Fill(1.);
869   }
870
871   // fill cut statistics and cut correlation histograms with information from the bitmap
872   if (f) return;
873
874   // Get the bitmap of the single cuts
875   SelectionBitMap(obj);
876
877   // Number of single cuts in this class
878   UInt_t ncuts = fBitmap->GetNbits();
879   for(UInt_t bit=0; bit<ncuts;bit++) {
880     if (!fBitmap->TestBitNumber(bit)) {
881         fhCutStatistics->Fill(bit+1);
882         for (UInt_t bit2=bit; bit2<ncuts;bit2++) {
883           if (!fBitmap->TestBitNumber(bit2)) 
884             fhCutCorrelation->Fill(bit+1,bit2+1);
885         }
886     }
887   }
888 }
889 //__________________________________________________________________________________
890 void AliCFTrackIsPrimaryCuts::SaveHistograms(const Char_t* dir) {
891   //
892   // saves the histograms in a directory (dir)
893   //
894   if(!fIsQAOn) return;
895
896   if (!dir)
897     dir = GetName();
898
899   gDirectory->mkdir(dir);
900   gDirectory->cd(dir);
901
902   gDirectory->mkdir("before_cuts");
903   gDirectory->mkdir("after_cuts");
904
905   fhCutStatistics->Write();
906   fhCutCorrelation->Write();
907
908   for (Int_t j=0; j<kNStepQA; j++) {
909     if (j==0)
910       gDirectory->cd("before_cuts");
911     else
912       gDirectory->cd("after_cuts");
913
914     fhDcaXYvsDcaZ[j]    ->Write();
915
916     for(Int_t i=0; i<kNHist; i++) fhQA[i][j]->Write();
917
918     gDirectory->cd("../");
919   }
920   gDirectory->cd("../");
921 }
922 //__________________________________________________________________________________
923 void AliCFTrackIsPrimaryCuts::DrawHistograms()
924 {
925   //
926   // draws some histograms
927   //
928   if(!fIsQAOn) return;
929
930   // pad margins
931   Float_t right = 0.03;
932   Float_t left = 0.175;
933   Float_t top = 0.03;
934   Float_t bottom = 0.175;
935
936   TCanvas* canvas1 = new TCanvas("Track_QA_Primary_1", "Track QA Primary 1", 800, 500);
937   canvas1->Divide(2, 1);
938
939   canvas1->cd(1);
940   fhCutStatistics->SetStats(kFALSE);
941   fhCutStatistics->LabelsOption("v");
942   gPad->SetLeftMargin(left);
943   gPad->SetBottomMargin(0.25);
944   gPad->SetRightMargin(right);
945   gPad->SetTopMargin(0.1);
946   fhCutStatistics->Draw();
947
948   canvas1->cd(2);
949   fhCutCorrelation->SetStats(kFALSE);
950   fhCutCorrelation->LabelsOption("v");
951   gPad->SetLeftMargin(0.30);
952   gPad->SetRightMargin(bottom);
953   gPad->SetTopMargin(0.1);
954   gPad->SetBottomMargin(0.25);
955   fhCutCorrelation->Draw("COLZ");
956
957   canvas1->SaveAs(Form("%s.eps", canvas1->GetName()));
958   canvas1->SaveAs(Form("%s.ps", canvas1->GetName()));
959
960   // -----
961
962   TCanvas* canvas2 = new TCanvas("Track_QA_Primary_2", "Track QA Primary 2", 800, 800);
963   canvas2->Divide(2, 2);
964
965   canvas2->cd(1);
966   gPad->SetRightMargin(right);
967   gPad->SetLeftMargin(left);
968   gPad->SetTopMargin(top);
969   gPad->SetBottomMargin(bottom);
970   fhQA[kDcaXY][0]->SetStats(kFALSE);
971   fhQA[kDcaXY][0]->Draw();
972   fhQA[kDcaXY][1]->Draw("same");
973
974   canvas2->cd(2);
975   gPad->SetRightMargin(right);
976   gPad->SetLeftMargin(left);
977   gPad->SetTopMargin(top);
978   gPad->SetBottomMargin(bottom);
979   fhQA[kDcaZ][0]->SetStats(kFALSE);
980   fhQA[kDcaZ][0]->Draw();
981   fhQA[kDcaZ][1]->Draw("same");
982
983   canvas2->cd(3);
984 //   fhDXYvsDZ[0]->SetStats(kFALSE);
985   gPad->SetLogz();
986   gPad->SetLeftMargin(bottom);
987   gPad->SetTopMargin(0.1);
988   gPad->SetBottomMargin(bottom);
989   gPad->SetRightMargin(0.2);
990   fhDcaXYvsDcaZ[0]->Draw("COLZ");
991
992   canvas2->cd(4);
993 //   fhDXYvsDZ[1]->SetStats(kFALSE);
994   gPad->SetLogz();
995   gPad->SetLeftMargin(bottom);
996   gPad->SetTopMargin(0.1);
997   gPad->SetBottomMargin(bottom);
998   gPad->SetRightMargin(0.2);
999   fhDcaXYvsDcaZ[1]->Draw("COLZ");
1000
1001   canvas2->SaveAs(Form("%s.eps", canvas2->GetName()));
1002   canvas2->SaveAs(Form("%s.ps", canvas2->GetName()));
1003
1004   // -----
1005
1006   TCanvas* canvas3 = new TCanvas("Track_QA_Primary_3", "Track QA Primary 3", 800, 400);
1007   canvas3->Divide(2, 1);
1008
1009   canvas3->cd(1);
1010   gPad->SetRightMargin(right);
1011   gPad->SetLeftMargin(left);
1012   gPad->SetTopMargin(top);
1013   gPad->SetBottomMargin(bottom);
1014   fhQA[kDcaXYnorm][0]->SetStats(kFALSE);
1015   fhQA[kDcaXYnorm][0]->Draw();
1016   fhQA[kDcaXYnorm][1]->Draw("same");
1017
1018   canvas3->cd(2);
1019   gPad->SetRightMargin(right);
1020   gPad->SetLeftMargin(left);
1021   gPad->SetTopMargin(top);
1022   gPad->SetBottomMargin(bottom);
1023   fhQA[kDcaZnorm][0]->SetStats(kFALSE);
1024   fhQA[kDcaZnorm][0]->Draw();
1025   fhQA[kDcaZnorm][1]->Draw("same");
1026
1027   canvas3->SaveAs(Form("%s.eps", canvas3->GetName()));
1028   canvas3->SaveAs(Form("%s.ps", canvas3->GetName()));
1029
1030   // -----
1031
1032   TCanvas* canvas4 = new TCanvas("Track_QA_Primary_4", "Track QA Primary 4", 1200, 800);
1033   canvas4->Divide(3, 2);
1034
1035   canvas4->cd(1);
1036   gPad->SetRightMargin(right);
1037   gPad->SetLeftMargin(left);
1038   gPad->SetTopMargin(top);
1039   gPad->SetBottomMargin(bottom);
1040   fhQA[kSigmaDcaXY][0]->SetStats(kFALSE);
1041   fhQA[kSigmaDcaXY][0]->Draw();
1042   fhQA[kSigmaDcaXY][1]->Draw("same");
1043
1044   canvas4->cd(2);
1045   gPad->SetRightMargin(right);
1046   gPad->SetLeftMargin(left);
1047   gPad->SetTopMargin(top);
1048   gPad->SetBottomMargin(bottom);
1049   fhQA[kSigmaDcaZ][0]->SetStats(kFALSE);
1050   fhQA[kSigmaDcaZ][0]->Draw();
1051   fhQA[kSigmaDcaZ][1]->Draw("same");
1052
1053   canvas4->cd(4);
1054   gPad->SetRightMargin(right);
1055   gPad->SetLeftMargin(left);
1056   gPad->SetTopMargin(top);
1057   gPad->SetBottomMargin(bottom);
1058   fhQA[kCutNSigmaToVertex][0]->SetStats(kFALSE);
1059   fhQA[kCutNSigmaToVertex][0]->Draw();
1060   fhQA[kCutNSigmaToVertex][1]->Draw("same");
1061
1062   canvas4->cd(5);
1063   gPad->SetRightMargin(right);
1064   gPad->SetLeftMargin(left);
1065   gPad->SetTopMargin(top);
1066   gPad->SetBottomMargin(bottom);
1067   fhQA[kCutRequireSigmaToVertex][0]->SetStats(kFALSE);
1068   fhQA[kCutRequireSigmaToVertex][0]->Draw();
1069   fhQA[kCutRequireSigmaToVertex][1]->Draw("same");
1070
1071   canvas4->cd(6);
1072   gPad->SetRightMargin(right);
1073   gPad->SetLeftMargin(left);
1074   gPad->SetTopMargin(top);
1075   gPad->SetBottomMargin(bottom);
1076   fhQA[kCutAcceptKinkDaughters][0]->SetStats(kFALSE);
1077   fhQA[kCutAcceptKinkDaughters][0]->Draw();
1078   fhQA[kCutAcceptKinkDaughters][1]->Draw("same");
1079
1080   canvas4->SaveAs(Form("%s.eps", canvas4->GetName()));
1081   canvas4->SaveAs(Form("%s.ps", canvas4->GetName()));
1082 }
1083 //__________________________________________________________________________________
1084 void AliCFTrackIsPrimaryCuts::AddQAHistograms(TList *qaList) {
1085   //
1086   // saves the histograms in a TList
1087   //
1088   DefineHistograms();
1089
1090   qaList->Add(fhCutStatistics);
1091   qaList->Add(fhCutCorrelation);
1092
1093   for (Int_t j=0; j<kNStepQA; j++) {
1094     qaList->Add(fhDcaXYvsDcaZ[j]);
1095     for(Int_t i=0; i<kNHist; i++)
1096         qaList->Add(fhQA[i][j]);
1097   }
1098 }