1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
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
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
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
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.
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.
46 // author: I. Kraus (Ingrid.Kraus@cern.ch)
48 // AliESDtrackCuts writte by Jan Fiete Grosse-Oetringhaus and
49 // AliRsnDaughterCut class written by A. Pulvirenti.
52 #include <TDirectory.h>
56 #include <AliESDtrack.h>
57 #include <AliESDEvent.h>
59 #include "AliCFTrackIsPrimaryCuts.h"
61 ClassImp(AliCFTrackIsPrimaryCuts)
63 //__________________________________________________________________________________
64 AliCFTrackIsPrimaryCuts::AliCFTrackIsPrimaryCuts() :
75 fNSigmaToVertexMin(0),
76 fNSigmaToVertexMax(0),
79 fRequireSigmaToVertex(0),
80 fAODType(AliAODTrack::kUndef),
81 fAcceptKinkDaughters(0),
86 fhNBinsRequireSigma(0),
95 fhBinLimRequireSigma(0x0),
96 fhBinLimAcceptKink(0x0),
99 fhBinLimDcaXYnorm(0x0),
100 fhBinLimDcaZnorm(0x0),
101 fhBinLimSigmaDcaXY(0x0),
102 fhBinLimSigmaDcaZ(0x0)
105 // Default constructor
109 //__________________________________________________________________________________
110 AliCFTrackIsPrimaryCuts::AliCFTrackIsPrimaryCuts(Char_t* name, Char_t* title) :
111 AliCFCutBase(name,title),
115 fMinDCAToVertexXY(0),
117 fMaxDCAToVertexXY(0),
121 fNSigmaToVertexMin(0),
122 fNSigmaToVertexMax(0),
125 fRequireSigmaToVertex(0),
126 fAODType(AliAODTrack::kUndef),
127 fAcceptKinkDaughters(0),
132 fhNBinsRequireSigma(0),
133 fhNBinsAcceptKink(0),
138 fhNBinsSigmaDcaXY(0),
141 fhBinLimRequireSigma(0x0),
142 fhBinLimAcceptKink(0x0),
145 fhBinLimDcaXYnorm(0x0),
146 fhBinLimDcaZnorm(0x0),
147 fhBinLimSigmaDcaXY(0x0),
148 fhBinLimSigmaDcaZ(0x0)
155 //__________________________________________________________________________________
156 AliCFTrackIsPrimaryCuts::AliCFTrackIsPrimaryCuts(const AliCFTrackIsPrimaryCuts& c) :
159 fUseSPDvertex(c.fUseSPDvertex),
160 fUseTPCvertex(c.fUseTPCvertex),
161 fMinDCAToVertexXY(c.fMinDCAToVertexXY),
162 fMinDCAToVertexZ(c.fMinDCAToVertexZ),
163 fMaxDCAToVertexXY(c.fMaxDCAToVertexXY),
164 fMaxDCAToVertexZ(c.fMaxDCAToVertexZ),
165 fDCAToVertex2D(c.fDCAToVertex2D),
166 fAbsDCAToVertex(c.fAbsDCAToVertex),
167 fNSigmaToVertexMin(c.fNSigmaToVertexMin),
168 fNSigmaToVertexMax(c.fNSigmaToVertexMax),
169 fSigmaDCAxy(c.fSigmaDCAxy),
170 fSigmaDCAz(c.fSigmaDCAz),
171 fRequireSigmaToVertex(c.fRequireSigmaToVertex),
172 fAODType(c.fAODType),
173 fAcceptKinkDaughters(c.fAcceptKinkDaughters),
174 fhCutStatistics(c.fhCutStatistics),
175 fhCutCorrelation(c.fhCutCorrelation),
177 fhNBinsNSigma(c.fhNBinsNSigma),
178 fhNBinsRequireSigma(c.fhNBinsRequireSigma),
179 fhNBinsAcceptKink(c.fhNBinsAcceptKink),
180 fhNBinsDcaXY(c.fhNBinsDcaXY),
181 fhNBinsDcaZ(c.fhNBinsDcaZ),
182 fhNBinsDcaXYnorm(c.fhNBinsDcaXYnorm),
183 fhNBinsDcaZnorm(c.fhNBinsDcaZnorm),
184 fhNBinsSigmaDcaXY(c.fhNBinsSigmaDcaXY),
185 fhNBinsSigmaDcaZ(c.fhNBinsSigmaDcaZ),
186 fhBinLimNSigma(c.fhBinLimNSigma),
187 fhBinLimRequireSigma(c.fhBinLimRequireSigma),
188 fhBinLimAcceptKink(c.fhBinLimAcceptKink),
189 fhBinLimDcaXY(c.fhBinLimDcaXY),
190 fhBinLimDcaZ(c.fhBinLimDcaZ),
191 fhBinLimDcaXYnorm(c.fhBinLimDcaXYnorm),
192 fhBinLimDcaZnorm(c.fhBinLimDcaZnorm),
193 fhBinLimSigmaDcaXY(c.fhBinLimSigmaDcaXY),
194 fhBinLimSigmaDcaZ(c.fhBinLimSigmaDcaZ)
199 ((AliCFTrackIsPrimaryCuts &) c).Copy(*this);
201 //__________________________________________________________________________________
202 AliCFTrackIsPrimaryCuts& AliCFTrackIsPrimaryCuts::operator=(const AliCFTrackIsPrimaryCuts& c)
205 // Assignment operator
208 AliCFCutBase::operator=(c) ;
210 fUseSPDvertex = c.fUseSPDvertex;
211 fUseTPCvertex = c.fUseTPCvertex;
212 fMinDCAToVertexXY = c.fMinDCAToVertexXY;
213 fMinDCAToVertexZ = c.fMinDCAToVertexZ;
214 fMaxDCAToVertexXY = c.fMaxDCAToVertexXY;
215 fMaxDCAToVertexZ = c.fMaxDCAToVertexZ;
216 fDCAToVertex2D = c.fDCAToVertex2D;
217 fAbsDCAToVertex = c.fAbsDCAToVertex;
218 fNSigmaToVertexMin = c.fNSigmaToVertexMin ;
219 fNSigmaToVertexMax = c.fNSigmaToVertexMax ;
220 fSigmaDCAxy = c.fSigmaDCAxy ;
221 fSigmaDCAz = c.fSigmaDCAz ;
222 fRequireSigmaToVertex = c.fRequireSigmaToVertex ;
223 fAODType = c.fAODType ;
224 fAcceptKinkDaughters = c.fAcceptKinkDaughters ;
225 fhCutStatistics = c.fhCutStatistics ;
226 fhCutCorrelation = c.fhCutCorrelation ;
228 fhNBinsNSigma = c.fhNBinsNSigma;
229 fhNBinsRequireSigma = c.fhNBinsRequireSigma;
230 fhNBinsAcceptKink = c.fhNBinsAcceptKink;
231 fhNBinsDcaXY = c.fhNBinsDcaXY;
232 fhNBinsDcaZ = c.fhNBinsDcaZ;
233 fhNBinsDcaXYnorm = c.fhNBinsDcaXYnorm;
234 fhNBinsDcaZnorm = c.fhNBinsDcaZnorm;
235 fhNBinsSigmaDcaXY = c.fhNBinsSigmaDcaXY;
236 fhNBinsSigmaDcaZ = c.fhNBinsSigmaDcaZ;
237 fhBinLimNSigma = c.fhBinLimNSigma;
238 fhBinLimRequireSigma = c.fhBinLimRequireSigma;
239 fhBinLimAcceptKink = c.fhBinLimAcceptKink;
240 fhBinLimDcaXY = c.fhBinLimDcaXY;
241 fhBinLimDcaZ = c.fhBinLimDcaZ;
242 fhBinLimDcaXYnorm = c.fhBinLimDcaXYnorm;
243 fhBinLimDcaZnorm = c.fhBinLimDcaZnorm;
244 fhBinLimSigmaDcaXY = c.fhBinLimSigmaDcaXY;
245 fhBinLimSigmaDcaZ = c.fhBinLimSigmaDcaZ;
247 for (Int_t j=0; j<6; j++) fDCA[j] = c.fDCA[j];
248 for (Int_t j=0; j<c.kNStepQA; j++){
249 if(c.fhDcaXYvsDcaZ[j]) fhDcaXYvsDcaZ[j] = (TH2F*)c.fhDcaXYvsDcaZ[j]->Clone();
250 for (Int_t i=0; i<c.kNHist; i++){
251 if(c.fhQA[i][j]) fhQA[i][j] = (TH1F*)c.fhQA[i][j]->Clone();
254 ((AliCFTrackIsPrimaryCuts &) c).Copy(*this);
258 //__________________________________________________________________________________
259 AliCFTrackIsPrimaryCuts::~AliCFTrackIsPrimaryCuts()
264 if (fhCutStatistics) delete fhCutStatistics;
265 if (fhCutCorrelation) delete fhCutCorrelation;
267 for (Int_t j=0; j<kNStepQA; j++){
268 if(fhDcaXYvsDcaZ[j]) delete fhDcaXYvsDcaZ[j];
269 for (Int_t i=0; i<kNHist; i++)
270 if(fhQA[i][j]) delete fhQA[i][j];
272 if(fESD) delete fESD;
273 if(fBitmap) delete fBitmap;
274 if(fhBinLimNSigma) delete fhBinLimNSigma;
275 if(fhBinLimRequireSigma) delete fhBinLimRequireSigma;
276 if(fhBinLimAcceptKink) delete fhBinLimAcceptKink;
277 if(fhBinLimDcaXY) delete fhBinLimDcaXY;
278 if(fhBinLimDcaZ) delete fhBinLimDcaZ;
279 if(fhBinLimDcaXYnorm) delete fhBinLimDcaXYnorm;
280 if(fhBinLimDcaZnorm) delete fhBinLimDcaZnorm;
281 if(fhBinLimSigmaDcaXY) delete fhBinLimSigmaDcaXY;
282 if(fhBinLimSigmaDcaZ) delete fhBinLimSigmaDcaZ;
284 //__________________________________________________________________________________
285 void AliCFTrackIsPrimaryCuts::Initialise()
288 // sets everything to zero
292 fMinDCAToVertexXY = 0;
293 fMinDCAToVertexZ = 0;
294 fMaxDCAToVertexXY = 0;
295 fMaxDCAToVertexZ = 0;
298 fNSigmaToVertexMin = 0;
299 fNSigmaToVertexMax = 0;
302 fRequireSigmaToVertex = 0;
303 fAcceptKinkDaughters = 0;
304 fAODType = AliAODTrack::kUndef;
306 SetMinDCAToVertexXY();
307 SetMinDCAToVertexZ();
308 SetMaxDCAToVertexXY();
309 SetMaxDCAToVertexZ();
312 SetMinNSigmaToVertex();
313 SetMaxNSigmaToVertex();
316 SetRequireSigmaToVertex();
317 SetAcceptKinkDaughters();
320 for (Int_t j=0; j<6; j++) fDCA[j] = 0.;
321 for (Int_t j=0; j<kNStepQA; j++) {
322 fhDcaXYvsDcaZ[j] = 0x0;
323 for (Int_t i=0; i<kNHist; i++)
327 fhCutCorrelation = 0;
328 fBitmap=new TBits(0);
330 //set default bining for QA histograms
331 SetHistogramBins(kCutNSigmaToVertex,100,0.,10.);
332 SetHistogramBins(kCutRequireSigmaToVertex,5,-0.75,1.75);
333 SetHistogramBins(kCutAcceptKinkDaughters,5,-0.75,1.75);
334 SetHistogramBins(kDcaXY,500,-10.,10.);
335 SetHistogramBins(kDcaZ,500,-10.,10.);
336 SetHistogramBins(kDcaXYnorm,500,-10.,10.);
337 SetHistogramBins(kDcaZnorm,500,-10.,10.);
338 SetHistogramBins(kSigmaDcaXY,500,-0.1,0.9);
339 SetHistogramBins(kSigmaDcaZ,500,-0.1,0.9);
341 //__________________________________________________________________________________
342 void AliCFTrackIsPrimaryCuts::Copy(TObject &c) const
347 AliCFTrackIsPrimaryCuts& target = (AliCFTrackIsPrimaryCuts &) c;
351 if (fhCutStatistics) target.fhCutStatistics = (TH1F*) fhCutStatistics->Clone();
352 if (fhCutCorrelation) target.fhCutCorrelation = (TH2F*) fhCutCorrelation->Clone();
354 for (Int_t j=0; j<6; j++) target.fDCA[j] = fDCA[j];
355 for (Int_t j=0; j<kNStepQA; j++){
356 if(fhDcaXYvsDcaZ[j]) target.fhDcaXYvsDcaZ[j] = (TH2F*)fhDcaXYvsDcaZ[j]->Clone();
357 for (Int_t i=0; i<kNHist; i++)
358 if(fhQA[i][j]) target.fhQA[i][j] = (TH1F*)fhQA[i][j]->Clone();
362 //__________________________________________________________________________________
363 void AliCFTrackIsPrimaryCuts::SetEvtInfo(TObject* esd) {
365 // Sets pointer to esd event information (AliESDEvent)
368 AliError("Pointer to AliESDEvent !");
371 TString className(esd->ClassName());
372 if (className.CompareTo("AliESDEvent") != 0) {
373 AliError("argument must point to an AliESDEvent !");
376 fESD = (AliESDEvent*) esd;
378 //__________________________________________________________________________________
379 void AliCFTrackIsPrimaryCuts::UseSPDvertex(Bool_t b) {
381 if(fUseTPCvertex && fUseSPDvertex) {
382 fUseSPDvertex = kFALSE;
383 AliError("SPD and TPC vertex chosen. TPC vertex is preferred.");
386 //__________________________________________________________________________________
387 void AliCFTrackIsPrimaryCuts::UseTPCvertex(Bool_t b) {
389 if(fUseTPCvertex) fUseSPDvertex = kFALSE;
391 //__________________________________________________________________________________
392 void AliCFTrackIsPrimaryCuts::GetDCA(AliESDtrack* esdTrack)
394 if (!esdTrack) return;
396 Float_t b[2] = {0.,0.};
397 Float_t bCov[3] = {0.,0.,0.};
398 if(!fUseSPDvertex && !fUseTPCvertex) esdTrack->GetImpactParameters(b,bCov);
399 if( fUseTPCvertex) esdTrack->GetImpactParametersTPC(b,bCov);
403 const AliESDVertex *vtx = fESD->GetVertex();
404 const Double_t Bz = fESD->GetMagneticField();
405 AliExternalTrackParam *cParam = 0;
406 Bool_t success = esdTrack->RelateToVertex(vtx, Bz, kVeryBig, cParam);
407 if (success) esdTrack->GetImpactParameters(b,bCov);
410 if (bCov[0]<=0 || bCov[2]<=0) {
411 bCov[0]=0; bCov[2]=0;
413 fDCA[0] = b[0]; // impact parameter xy
414 fDCA[1] = b[1]; // impact parameter z
415 fDCA[2] = TMath::Sqrt(bCov[0]); // resolution xy
416 fDCA[3] = TMath::Sqrt(bCov[2]); // resolution z
418 if (!fAbsDCAToVertex) {
419 if (fDCA[2] > 0) fDCA[0] = fDCA[0]/fDCA[2]; // normalised impact parameter xy
420 if (fDCA[3] > 0) fDCA[1] = fDCA[1]/fDCA[3]; // normalised impact parameter z
424 if(!fUseSPDvertex && !fUseTPCvertex)
425 fDCA[5] = AliESDtrackCuts::GetSigmaToVertex(esdTrack);
429 if (fDCA[2]==0 || fDCA[3]==0)
432 Float_t d = TMath::Sqrt(TMath::Power(b[0]/fDCA[2],2) + TMath::Power(b[1]/fDCA[3],2));
433 if (TMath::Exp(-d * d / 2) < 1e-15)
435 fDCA[5] = TMath::ErfInverse(1 - TMath::Exp(-d * d / 2)) * TMath::Sqrt(2);
439 //__________________________________________________________________________________
440 void AliCFTrackIsPrimaryCuts::SelectionBitMap(TObject* obj)
443 // test if the track passes the single cuts
444 // and store the information in a bitmap
447 // bitmap stores the decision of each single cut
448 for(Int_t i=0; i<kNCuts; i++)fBitmap->SetBitNumber(i,kFALSE);
450 // check TObject and cast into ESDtrack
452 if (!obj->InheritsFrom("AliVParticle")) {
453 AliError("object must derived from AliVParticle !");
457 Bool_t isESDTrack = strcmp(obj->ClassName(),"AliESDtrack") == 0 ? kTRUE : kFALSE ;
458 Bool_t isAODTrack = strcmp(obj->ClassName(),"AliAODTrack") == 0 ? kTRUE : kFALSE ;
460 AliESDtrack * esdTrack = 0x0 ;
461 AliAODTrack * aodTrack = 0x0 ;
462 if (isESDTrack) esdTrack = dynamic_cast<AliESDtrack*>(obj);
463 if (isAODTrack) aodTrack = dynamic_cast<AliAODTrack*>(obj);
465 // get the track to vertex parameter for ESD track
466 if (isESDTrack) GetDCA(esdTrack);
468 Float_t bxy = 0, bz = 0;
470 bxy = TMath::Abs(fDCA[0]);
471 bz = TMath::Abs(fDCA[1]);
473 Float_t b2Dmin = 0, b2Dmax = 0;
475 if (fMinDCAToVertexXY>0 && fMinDCAToVertexZ>0)
476 b2Dmin = fDCA[0]*fDCA[0]/fMinDCAToVertexXY/fMinDCAToVertexXY + fDCA[1]*fDCA[1]/fMinDCAToVertexZ/fMinDCAToVertexZ;
477 if (fMaxDCAToVertexXY>0 && fMaxDCAToVertexZ>0)
478 b2Dmax = fDCA[0]*fDCA[0]/fMaxDCAToVertexXY/fMaxDCAToVertexXY + fDCA[1]*fDCA[1]/fMaxDCAToVertexZ/fMaxDCAToVertexZ;
485 if (fDCAToVertex2D || (!fDCAToVertex2D && bxy >= fMinDCAToVertexXY && bxy <= fMaxDCAToVertexXY)) {
486 fBitmap->SetBitNumber(iCutBit,kTRUE);
489 else fBitmap->SetBitNumber(iCutBit,kTRUE);
494 if (fDCAToVertex2D || (!fDCAToVertex2D && bz >= fMinDCAToVertexZ && bz <= fMaxDCAToVertexZ)) {
495 fBitmap->SetBitNumber(iCutBit,kTRUE);
498 else fBitmap->SetBitNumber(iCutBit,kTRUE);
503 if (!fDCAToVertex2D || (fDCAToVertex2D && b2Dmin > 1 && b2Dmax < 1)) {
504 fBitmap->SetBitNumber(iCutBit,kTRUE);
507 else fBitmap->SetBitNumber(iCutBit,kTRUE);
512 if (fDCA[5] >= fNSigmaToVertexMin && fDCA[5] <= fNSigmaToVertexMax) {
513 fBitmap->SetBitNumber(iCutBit,kTRUE);
516 else fBitmap->SetBitNumber(iCutBit,kTRUE);
521 if (fDCA[2] < fSigmaDCAxy) {
522 fBitmap->SetBitNumber(iCutBit,kTRUE);
525 else fBitmap->SetBitNumber(iCutBit,kTRUE);
530 if (fDCA[3] < fSigmaDCAz) {
531 fBitmap->SetBitNumber(iCutBit,kTRUE);
534 else fBitmap->SetBitNumber(iCutBit,kTRUE);
539 if (!fRequireSigmaToVertex || (fDCA[5]>=0 && fRequireSigmaToVertex)) {
540 fBitmap->SetBitNumber(iCutBit,kTRUE);
543 else fBitmap->SetBitNumber(iCutBit,kTRUE);
548 if (fAcceptKinkDaughters || (!fAcceptKinkDaughters && esdTrack->GetKinkIndex(0)<=0)) {
549 fBitmap->SetBitNumber(iCutBit,kTRUE);
552 else fBitmap->SetBitNumber(iCutBit,kTRUE);
557 if (fAODType==AliAODTrack::kUndef || fAODType == aodTrack->GetType()) {
558 fBitmap->SetBitNumber(iCutBit,kTRUE);
561 else fBitmap->SetBitNumber(iCutBit,kTRUE);
565 //__________________________________________________________________________________
566 Bool_t AliCFTrackIsPrimaryCuts::IsSelected(TObject* obj) {
568 // loops over decisions of single cuts and returns if the track is accepted
570 SelectionBitMap(obj);
572 if (fIsQAOn) FillHistograms(obj,0);
573 Bool_t isSelected = kTRUE;
575 for (UInt_t icut=0; icut<fBitmap->GetNbits();icut++) {
576 if(!fBitmap->TestBitNumber(icut)) {
581 if (!isSelected) return kFALSE ;
582 if (fIsQAOn) FillHistograms(obj,1);
585 //__________________________________________________________________________________
586 void AliCFTrackIsPrimaryCuts::SetHistogramBins(Int_t index, Int_t nbins, Double_t *bins)
594 case kCutNSigmaToVertex:
595 fhNBinsNSigma=nbins+1;
596 fhBinLimNSigma=new Double_t[nbins+1];
597 for(Int_t i=0;i<nbins+1;i++)fhBinLimNSigma[i]=bins[i];
600 case kCutRequireSigmaToVertex:
601 fhNBinsRequireSigma=nbins+1;
602 fhBinLimRequireSigma=new Double_t[nbins+1];
603 for(Int_t i=0;i<nbins+1;i++)fhBinLimRequireSigma[i]=bins[i];
606 case kCutAcceptKinkDaughters:
607 fhNBinsAcceptKink=nbins+1;
608 fhBinLimAcceptKink=new Double_t[nbins+1];
609 for(Int_t i=0;i<nbins+1;i++)fhBinLimAcceptKink[i]=bins[i];
613 fhNBinsDcaXY=nbins+1;
614 fhBinLimDcaXY=new Double_t[nbins+1];
615 for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaXY[i]=bins[i];
620 fhBinLimDcaZ=new Double_t[nbins+1];
621 for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaZ[i]=bins[i];
625 fhNBinsDcaXYnorm=nbins+1;
626 fhBinLimDcaXYnorm=new Double_t[nbins+1];
627 for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaXYnorm[i]=bins[i];
631 fhNBinsDcaZnorm=nbins+1;
632 fhBinLimDcaZnorm=new Double_t[nbins+1];
633 for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaZnorm[i]=bins[i];
637 fhNBinsSigmaDcaXY=nbins+1;
638 fhBinLimSigmaDcaXY=new Double_t[nbins+1];
639 for(Int_t i=0;i<nbins+1;i++)fhBinLimSigmaDcaXY[i]=bins[i];
643 fhNBinsSigmaDcaZ=nbins+1;
644 fhBinLimSigmaDcaZ=new Double_t[nbins+1];
645 for(Int_t i=0;i<nbins+1;i++)fhBinLimSigmaDcaZ[i]=bins[i];
649 //__________________________________________________________________________________
650 void AliCFTrackIsPrimaryCuts::SetHistogramBins(Int_t index, Int_t nbins, Double_t xmin, Double_t xmax)
656 case kCutNSigmaToVertex:
657 fhNBinsNSigma=nbins+1;
658 fhBinLimNSigma=new Double_t[nbins+1];
659 for(Int_t i=0;i<nbins+1;i++)fhBinLimNSigma[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
662 case kCutRequireSigmaToVertex:
663 fhNBinsRequireSigma=nbins+1;
664 fhBinLimRequireSigma=new Double_t[nbins+1];
665 for(Int_t i=0;i<nbins+1;i++)fhBinLimRequireSigma[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
668 case kCutAcceptKinkDaughters:
669 fhNBinsAcceptKink=nbins+1;
670 fhBinLimAcceptKink=new Double_t[nbins+1];
671 for(Int_t i=0;i<nbins+1;i++)fhBinLimAcceptKink[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
675 fhNBinsDcaXY=nbins+1;
676 fhBinLimDcaXY=new Double_t[nbins+1];
677 for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaXY[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
682 fhBinLimDcaZ=new Double_t[nbins+1];
683 for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaZ[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
687 fhNBinsDcaXYnorm=nbins+1;
688 fhBinLimDcaXYnorm=new Double_t[nbins+1];
689 for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaXYnorm[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
693 fhNBinsDcaZnorm=nbins+1;
694 fhBinLimDcaZnorm=new Double_t[nbins+1];
695 for(Int_t i=0;i<nbins+1;i++)fhBinLimDcaZnorm[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
699 fhNBinsSigmaDcaXY=nbins+1;
700 fhBinLimSigmaDcaXY=new Double_t[nbins+1];
701 for(Int_t i=0;i<nbins+1;i++)fhBinLimSigmaDcaXY[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
705 fhNBinsSigmaDcaZ=nbins+1;
706 fhBinLimSigmaDcaZ=new Double_t[nbins+1];
707 for(Int_t i=0;i<nbins+1;i++)fhBinLimSigmaDcaZ[i]=xmin+i*(xmax-xmin)/Double_t(nbins);
711 //__________________________________________________________________________________
712 void AliCFTrackIsPrimaryCuts::DefineHistograms() {
714 // histograms for cut variables, cut statistics and cut correlations
718 // book cut statistics and cut correlation histograms
719 fhCutStatistics = new TH1F(Form("%s_cut_statistics",GetName()), Form("%s cut statistics",GetName()), kNCuts,0.5,kNCuts+0.5);
720 fhCutStatistics->SetLineWidth(2);
721 fhCutStatistics->GetXaxis()->SetBinLabel(1,"dca xy");
722 fhCutStatistics->GetXaxis()->SetBinLabel(2,"dca z");
723 fhCutStatistics->GetXaxis()->SetBinLabel(3,"dca ellipse");
724 fhCutStatistics->GetXaxis()->SetBinLabel(4,"n dca");
725 fhCutStatistics->GetXaxis()->SetBinLabel(5,"sigma dca xy");
726 fhCutStatistics->GetXaxis()->SetBinLabel(6,"sigma dca z");
727 fhCutStatistics->GetXaxis()->SetBinLabel(7,"require dca");
728 fhCutStatistics->GetXaxis()->SetBinLabel(8,"kink daughter");
729 fhCutStatistics->GetXaxis()->SetBinLabel(9,"AOD type");
731 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);
732 fhCutCorrelation->SetLineWidth(2);
733 fhCutCorrelation->GetXaxis()->SetBinLabel(1,"dca xy");
734 fhCutCorrelation->GetXaxis()->SetBinLabel(2,"dca z");
735 fhCutCorrelation->GetXaxis()->SetBinLabel(3,"dca ellipse");
736 fhCutCorrelation->GetXaxis()->SetBinLabel(4,"n dca");
737 fhCutCorrelation->GetXaxis()->SetBinLabel(5,"sigma dca xy");
738 fhCutCorrelation->GetXaxis()->SetBinLabel(6,"sigma dca z");
739 fhCutCorrelation->GetXaxis()->SetBinLabel(7,"require dca");
740 fhCutCorrelation->GetXaxis()->SetBinLabel(8,"kink daughter");
741 fhCutCorrelation->GetXaxis()->SetBinLabel(9,"AOD type");
743 fhCutCorrelation->GetYaxis()->SetBinLabel(1,"dca xy");
744 fhCutCorrelation->GetYaxis()->SetBinLabel(2,"dca z");
745 fhCutCorrelation->GetYaxis()->SetBinLabel(3,"dca ellipse");
746 fhCutCorrelation->GetYaxis()->SetBinLabel(4,"n dca");
747 fhCutCorrelation->GetYaxis()->SetBinLabel(5,"sigma dca xy");
748 fhCutCorrelation->GetYaxis()->SetBinLabel(6,"sigma dca z");
749 fhCutCorrelation->GetYaxis()->SetBinLabel(7,"require dca");
750 fhCutCorrelation->GetYaxis()->SetBinLabel(8,"kink daughter");
751 fhCutCorrelation->GetYaxis()->SetBinLabel(9,"AOD type");
753 // book QA histograms
755 for (Int_t i=0; i<kNStepQA; i++) {
756 if (i==0) sprintf(str," ");
757 else sprintf(str,"_cut");
759 fhDcaXYvsDcaZ[i] = new TH2F(Form("%s_dcaXYvsDcaZ%s",GetName(),str),"",200,-10,10,200,-10,10);
760 fhQA[kCutNSigmaToVertex][i] = new TH1F(Form("%s_nSigmaToVertex%s",GetName(),str),"",fhNBinsNSigma-1,fhBinLimNSigma);
761 fhQA[kCutRequireSigmaToVertex][i] = new TH1F(Form("%s_requireSigmaToVertex%s",GetName(),str),"",fhNBinsRequireSigma-1,fhBinLimRequireSigma);
762 fhQA[kCutAcceptKinkDaughters][i] = new TH1F(Form("%s_acceptKinkDaughters%s",GetName(),str),"",fhNBinsAcceptKink-1,fhBinLimAcceptKink);
763 fhQA[kDcaXY][i] = new TH1F(Form("%s_dcaXY%s",GetName(),str),"",fhNBinsDcaXY-1,fhBinLimDcaXY);
764 fhQA[kDcaZ][i] = new TH1F(Form("%s_dcaZ%s",GetName(),str),"",fhNBinsDcaZ-1,fhBinLimDcaZ);
765 fhQA[kDcaXYnorm][i] = new TH1F(Form("%s_dcaXYnorm%s",GetName(),str),"",fhNBinsDcaXYnorm-1,fhBinLimDcaXYnorm);
766 fhQA[kDcaZnorm][i] = new TH1F(Form("%s_dcaZnorm%s",GetName(),str),"",fhNBinsDcaZnorm-1,fhBinLimDcaZnorm);
767 fhQA[kSigmaDcaXY][i] = new TH1F(Form("%s_sigmaDcaXY%s",GetName(),str),"",fhNBinsSigmaDcaXY-1,fhBinLimSigmaDcaXY);
768 fhQA[kSigmaDcaZ][i] = new TH1F(Form("%s_sigmaDcaZ%s",GetName(),str),"",fhNBinsSigmaDcaZ-1,fhBinLimSigmaDcaZ);
770 fhDcaXYvsDcaZ[i]->SetXTitle("impact par. d_{z}");
771 fhDcaXYvsDcaZ[i]->SetYTitle("impact par. d_{xy}");
773 fhQA[kCutNSigmaToVertex][i]->SetXTitle("n #sigma to vertex");
774 fhQA[kCutRequireSigmaToVertex][i]->SetXTitle("require #sigma to vertex");
775 fhQA[kCutAcceptKinkDaughters][i]->SetXTitle("accept kink daughters");
776 fhQA[kDcaXY][i]->SetXTitle("impact par. d_{xy}");
777 fhQA[kDcaZ][i]->SetXTitle("impact par. d_{z}");
778 fhQA[kDcaXYnorm][i]->SetXTitle("norm. impact par. d_{xy} / #sigma_{xy}");
779 fhQA[kDcaZnorm][i]->SetXTitle("norm. impact par. d_{z} / #sigma_{z}");
780 fhQA[kSigmaDcaXY][i]->SetXTitle("impact par. resolution #sigma_{xy}");
781 fhQA[kSigmaDcaZ][i]->SetXTitle("impact par. resolution #sigma_{z}");
783 for(Int_t i=0; i<kNHist; i++) fhQA[i][1]->SetLineColor(color);
785 //__________________________________________________________________________________
786 void AliCFTrackIsPrimaryCuts::FillHistograms(TObject* obj, Bool_t f)
789 // fill the QA histograms
794 Bool_t isESDTrack = strcmp(obj->ClassName(),"AliESDtrack") == 0 ? kTRUE : kFALSE ;
795 Bool_t isAODTrack = strcmp(obj->ClassName(),"AliAODTrack") == 0 ? kTRUE : kFALSE ;
797 AliESDtrack * esdTrack = 0x0 ;
798 AliAODTrack * aodTrack = 0x0 ;
799 if (isESDTrack) esdTrack = dynamic_cast<AliESDtrack*>(obj);
800 if (isAODTrack) aodTrack = dynamic_cast<AliAODTrack*>(obj);
802 // f = 0: fill histograms before cuts
803 // f = 1: fill histograms after cuts
805 // get the track to vertex parameter for ESD track
808 fhQA[kDcaZ][f]->Fill(fDCA[1]);
809 fhQA[kDcaXY][f]->Fill(fDCA[0]);
810 fhDcaXYvsDcaZ[f]->Fill(fDCA[1],fDCA[0]);
811 fhQA[kSigmaDcaXY][f]->Fill(fDCA[2]);
812 fhQA[kSigmaDcaZ][f]->Fill(fDCA[3]);
813 // // // // // // // delete histograms
814 fhQA[kDcaZnorm][f]->Fill(fDCA[1]);
815 fhQA[kDcaXYnorm][f]->Fill(fDCA[0]);
817 fhQA[kCutNSigmaToVertex][f]->Fill(fDCA[5]);
818 if (fDCA[5]<0 && fRequireSigmaToVertex) fhQA[kCutRequireSigmaToVertex][f]->Fill(0.);
819 if (!(fDCA[5]<0 && fRequireSigmaToVertex)) fhQA[kCutRequireSigmaToVertex][f]->Fill(1.);
821 if (!fAcceptKinkDaughters && esdTrack->GetKinkIndex(0)>0) fhQA[kCutAcceptKinkDaughters][f]->Fill(0.);
822 if (!(!fAcceptKinkDaughters && esdTrack->GetKinkIndex(0)>0)) fhQA[kCutAcceptKinkDaughters][f]->Fill(0.);
825 // fill cut statistics and cut correlation histograms with information from the bitmap
828 // Get the bitmap of the single cuts
830 SelectionBitMap(obj);
832 // Number of single cuts in this class
833 UInt_t ncuts = fBitmap->GetNbits();
834 for(UInt_t bit=0; bit<ncuts;bit++) {
835 if (!fBitmap->TestBitNumber(bit)) {
836 fhCutStatistics->Fill(bit+1);
837 for (UInt_t bit2=bit; bit2<ncuts;bit2++) {
838 if (!fBitmap->TestBitNumber(bit2))
839 fhCutCorrelation->Fill(bit+1,bit2+1);
844 //__________________________________________________________________________________
845 void AliCFTrackIsPrimaryCuts::SaveHistograms(const Char_t* dir) {
847 // saves the histograms in a directory (dir)
854 gDirectory->mkdir(dir);
857 gDirectory->mkdir("before_cuts");
858 gDirectory->mkdir("after_cuts");
860 fhCutStatistics->Write();
861 fhCutCorrelation->Write();
863 for (Int_t j=0; j<kNStepQA; j++) {
865 gDirectory->cd("before_cuts");
867 gDirectory->cd("after_cuts");
869 fhDcaXYvsDcaZ[j] ->Write();
871 for(Int_t i=0; i<kNHist; i++) fhQA[i][j]->Write();
873 gDirectory->cd("../");
875 gDirectory->cd("../");
877 //__________________________________________________________________________________
878 void AliCFTrackIsPrimaryCuts::DrawHistograms()
881 // draws some histograms
886 Float_t right = 0.03;
887 Float_t left = 0.175;
889 Float_t bottom = 0.175;
891 TCanvas* canvas1 = new TCanvas("Track_QA_Primary_1", "Track QA Primary 1", 800, 500);
892 canvas1->Divide(2, 1);
895 fhCutStatistics->SetStats(kFALSE);
896 fhCutStatistics->LabelsOption("v");
897 gPad->SetLeftMargin(left);
898 gPad->SetBottomMargin(0.25);
899 gPad->SetRightMargin(right);
900 gPad->SetTopMargin(0.1);
901 fhCutStatistics->Draw();
904 fhCutCorrelation->SetStats(kFALSE);
905 fhCutCorrelation->LabelsOption("v");
906 gPad->SetLeftMargin(0.30);
907 gPad->SetRightMargin(bottom);
908 gPad->SetTopMargin(0.1);
909 gPad->SetBottomMargin(0.25);
910 fhCutCorrelation->Draw("COLZ");
912 canvas1->SaveAs(Form("%s.eps", canvas1->GetName()));
913 canvas1->SaveAs(Form("%s.ps", canvas1->GetName()));
917 TCanvas* canvas2 = new TCanvas("Track_QA_Primary_2", "Track QA Primary 2", 800, 800);
918 canvas2->Divide(2, 2);
921 gPad->SetRightMargin(right);
922 gPad->SetLeftMargin(left);
923 gPad->SetTopMargin(top);
924 gPad->SetBottomMargin(bottom);
925 fhQA[kDcaXY][0]->SetStats(kFALSE);
926 fhQA[kDcaXY][0]->Draw();
927 fhQA[kDcaXY][1]->Draw("same");
930 gPad->SetRightMargin(right);
931 gPad->SetLeftMargin(left);
932 gPad->SetTopMargin(top);
933 gPad->SetBottomMargin(bottom);
934 fhQA[kDcaZ][0]->SetStats(kFALSE);
935 fhQA[kDcaZ][0]->Draw();
936 fhQA[kDcaZ][1]->Draw("same");
939 // fhDXYvsDZ[0]->SetStats(kFALSE);
941 gPad->SetLeftMargin(bottom);
942 gPad->SetTopMargin(0.1);
943 gPad->SetBottomMargin(bottom);
944 gPad->SetRightMargin(0.2);
945 fhDcaXYvsDcaZ[0]->Draw("COLZ");
948 // fhDXYvsDZ[1]->SetStats(kFALSE);
950 gPad->SetLeftMargin(bottom);
951 gPad->SetTopMargin(0.1);
952 gPad->SetBottomMargin(bottom);
953 gPad->SetRightMargin(0.2);
954 fhDcaXYvsDcaZ[1]->Draw("COLZ");
956 canvas2->SaveAs(Form("%s.eps", canvas2->GetName()));
957 canvas2->SaveAs(Form("%s.ps", canvas2->GetName()));
961 TCanvas* canvas3 = new TCanvas("Track_QA_Primary_3", "Track QA Primary 3", 800, 400);
962 canvas3->Divide(2, 1);
965 gPad->SetRightMargin(right);
966 gPad->SetLeftMargin(left);
967 gPad->SetTopMargin(top);
968 gPad->SetBottomMargin(bottom);
969 fhQA[kDcaXYnorm][0]->SetStats(kFALSE);
970 fhQA[kDcaXYnorm][0]->Draw();
971 fhQA[kDcaXYnorm][1]->Draw("same");
974 gPad->SetRightMargin(right);
975 gPad->SetLeftMargin(left);
976 gPad->SetTopMargin(top);
977 gPad->SetBottomMargin(bottom);
978 fhQA[kDcaZnorm][0]->SetStats(kFALSE);
979 fhQA[kDcaZnorm][0]->Draw();
980 fhQA[kDcaZnorm][1]->Draw("same");
982 canvas3->SaveAs(Form("%s.eps", canvas3->GetName()));
983 canvas3->SaveAs(Form("%s.ps", canvas3->GetName()));
987 TCanvas* canvas4 = new TCanvas("Track_QA_Primary_4", "Track QA Primary 4", 1200, 800);
988 canvas4->Divide(3, 2);
991 gPad->SetRightMargin(right);
992 gPad->SetLeftMargin(left);
993 gPad->SetTopMargin(top);
994 gPad->SetBottomMargin(bottom);
995 fhQA[kSigmaDcaXY][0]->SetStats(kFALSE);
996 fhQA[kSigmaDcaXY][0]->Draw();
997 fhQA[kSigmaDcaXY][1]->Draw("same");
1000 gPad->SetRightMargin(right);
1001 gPad->SetLeftMargin(left);
1002 gPad->SetTopMargin(top);
1003 gPad->SetBottomMargin(bottom);
1004 fhQA[kSigmaDcaZ][0]->SetStats(kFALSE);
1005 fhQA[kSigmaDcaZ][0]->Draw();
1006 fhQA[kSigmaDcaZ][1]->Draw("same");
1009 gPad->SetRightMargin(right);
1010 gPad->SetLeftMargin(left);
1011 gPad->SetTopMargin(top);
1012 gPad->SetBottomMargin(bottom);
1013 fhQA[kCutNSigmaToVertex][0]->SetStats(kFALSE);
1014 fhQA[kCutNSigmaToVertex][0]->Draw();
1015 fhQA[kCutNSigmaToVertex][1]->Draw("same");
1018 gPad->SetRightMargin(right);
1019 gPad->SetLeftMargin(left);
1020 gPad->SetTopMargin(top);
1021 gPad->SetBottomMargin(bottom);
1022 fhQA[kCutRequireSigmaToVertex][0]->SetStats(kFALSE);
1023 fhQA[kCutRequireSigmaToVertex][0]->Draw();
1024 fhQA[kCutRequireSigmaToVertex][1]->Draw("same");
1027 gPad->SetRightMargin(right);
1028 gPad->SetLeftMargin(left);
1029 gPad->SetTopMargin(top);
1030 gPad->SetBottomMargin(bottom);
1031 fhQA[kCutAcceptKinkDaughters][0]->SetStats(kFALSE);
1032 fhQA[kCutAcceptKinkDaughters][0]->Draw();
1033 fhQA[kCutAcceptKinkDaughters][1]->Draw("same");
1035 canvas4->SaveAs(Form("%s.eps", canvas4->GetName()));
1036 canvas4->SaveAs(Form("%s.ps", canvas4->GetName()));
1038 //__________________________________________________________________________________
1039 void AliCFTrackIsPrimaryCuts::AddQAHistograms(TList *qaList) {
1041 // saves the histograms in a TList
1045 qaList->Add(fhCutStatistics);
1046 qaList->Add(fhCutCorrelation);
1048 for (Int_t j=0; j<kNStepQA; j++) {
1049 qaList->Add(fhDcaXYvsDcaZ[j]);
1050 for(Int_t i=0; i<kNHist; i++)
1051 qaList->Add(fhQA[i][j]);