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 /* $Id: AliESDtrackCuts.cxx 24534 2008-03-16 22:22:11Z fca $ */
18 #include "AliESDtrackCuts.h"
20 #include <AliESDtrack.h>
21 #include <AliESDVertex.h>
22 #include <AliESDEvent.h>
27 #include <TDirectory.h>
29 //____________________________________________________________________
30 ClassImp(AliESDtrackCuts)
33 const Char_t* AliESDtrackCuts::fgkCutNames[kNCuts] = {
38 "#Chi^{2}/cluster TPC",
39 "#Chi^{2}/cluster ITS",
55 "trk-to-vtx max dca 2D absolute",
56 "trk-to-vtx max dca xy absolute",
57 "trk-to-vtx max dca z absolute",
58 "trk-to-vtx min dca 2D absolute",
59 "trk-to-vtx min dca xy absolute",
60 "trk-to-vtx min dca z absolute",
61 "SPD cluster requirement",
62 "SDD cluster requirement",
63 "SSD cluster requirement"
66 //____________________________________________________________________
67 AliESDtrackCuts::AliESDtrackCuts(const Char_t* name, const Char_t* title) : AliAnalysisCuts(name,title),
68 fCutMinNClusterTPC(0),
69 fCutMinNClusterITS(0),
70 fCutMaxChi2PerClusterTPC(0),
71 fCutMaxChi2PerClusterITS(0),
77 fCutAcceptKinkDaughters(0),
78 fCutRequireTPCRefit(0),
79 fCutRequireITSRefit(0),
80 fCutNsigmaToVertex(0),
81 fCutSigmaToVertexRequired(0),
82 fCutMaxDCAToVertexXY(0),
83 fCutMaxDCAToVertexZ(0),
84 fCutMinDCAToVertexXY(0),
85 fCutMinDCAToVertexZ(0),
112 //##############################################################################
113 // setting default cuts
114 SetMinNClustersTPC();
115 SetMinNClustersITS();
116 SetMaxChi2PerClusterTPC();
117 SetMaxChi2PerClusterITS();
118 SetMaxCovDiagonalElements();
119 SetRequireTPCRefit();
120 SetRequireITSRefit();
121 SetAcceptKinkDaughters();
122 SetMaxNsigmaToVertex();
123 SetMaxDCAToVertexXY();
124 SetMaxDCAToVertexZ();
126 SetMinDCAToVertexXY();
127 SetMinDCAToVertexZ();
135 SetClusterRequirementITS(kSPD);
136 SetClusterRequirementITS(kSDD);
137 SetClusterRequirementITS(kSSD);
142 //_____________________________________________________________________________
143 AliESDtrackCuts::AliESDtrackCuts(const AliESDtrackCuts &c) : AliAnalysisCuts(c),
144 fCutMinNClusterTPC(0),
145 fCutMinNClusterITS(0),
146 fCutMaxChi2PerClusterTPC(0),
147 fCutMaxChi2PerClusterITS(0),
153 fCutAcceptKinkDaughters(0),
154 fCutRequireTPCRefit(0),
155 fCutRequireITSRefit(0),
156 fCutNsigmaToVertex(0),
157 fCutSigmaToVertexRequired(0),
158 fCutMaxDCAToVertexXY(0),
159 fCutMaxDCAToVertexZ(0),
160 fCutMinDCAToVertexXY(0),
161 fCutMinDCAToVertexZ(0),
162 fCutDCAToVertex2D(0),
186 ((AliESDtrackCuts &) c).Copy(*this);
189 AliESDtrackCuts::~AliESDtrackCuts()
195 for (Int_t i=0; i<2; i++) {
197 if (fhNClustersITS[i])
198 delete fhNClustersITS[i];
199 if (fhNClustersTPC[i])
200 delete fhNClustersTPC[i];
201 if (fhChi2PerClusterITS[i])
202 delete fhChi2PerClusterITS[i];
203 if (fhChi2PerClusterTPC[i])
204 delete fhChi2PerClusterTPC[i];
225 if (fhDXYNormalized[i])
226 delete fhDXYNormalized[i];
227 if (fhDZNormalized[i])
228 delete fhDZNormalized[i];
229 if (fhDXYvsDZNormalized[i])
230 delete fhDXYvsDZNormalized[i];
231 if (fhNSigmaToVertex[i])
232 delete fhNSigmaToVertex[i];
240 delete ffDTheoretical;
243 delete fhCutStatistics;
244 if (fhCutCorrelation)
245 delete fhCutCorrelation;
248 void AliESDtrackCuts::Init()
251 // sets everything to zero
254 fCutMinNClusterTPC = 0;
255 fCutMinNClusterITS = 0;
257 fCutMaxChi2PerClusterTPC = 0;
258 fCutMaxChi2PerClusterITS = 0;
260 for (Int_t i = 0; i < 3; i++)
261 fCutClusterRequirementITS[i] = kOff;
269 fCutAcceptKinkDaughters = 0;
270 fCutRequireTPCRefit = 0;
271 fCutRequireITSRefit = 0;
273 fCutNsigmaToVertex = 0;
274 fCutSigmaToVertexRequired = 0;
275 fCutMaxDCAToVertexXY = 0;
276 fCutMaxDCAToVertexZ = 0;
277 fCutDCAToVertex2D = 0;
278 fCutMinDCAToVertexXY = 0;
279 fCutMinDCAToVertexZ = 0;
297 fHistogramsOn = kFALSE;
299 for (Int_t i=0; i<2; ++i)
301 fhNClustersITS[i] = 0;
302 fhNClustersTPC[i] = 0;
304 fhChi2PerClusterITS[i] = 0;
305 fhChi2PerClusterTPC[i] = 0;
318 fhDXYNormalized[i] = 0;
319 fhDZNormalized[i] = 0;
320 fhDXYvsDZNormalized[i] = 0;
321 fhNSigmaToVertex[i] = 0;
329 fhCutCorrelation = 0;
332 //_____________________________________________________________________________
333 AliESDtrackCuts &AliESDtrackCuts::operator=(const AliESDtrackCuts &c)
336 // Assignment operator
339 if (this != &c) ((AliESDtrackCuts &) c).Copy(*this);
343 //_____________________________________________________________________________
344 void AliESDtrackCuts::Copy(TObject &c) const
350 AliESDtrackCuts& target = (AliESDtrackCuts &) c;
354 target.fCutMinNClusterTPC = fCutMinNClusterTPC;
355 target.fCutMinNClusterITS = fCutMinNClusterITS;
357 target.fCutMaxChi2PerClusterTPC = fCutMaxChi2PerClusterTPC;
358 target.fCutMaxChi2PerClusterITS = fCutMaxChi2PerClusterITS;
360 for (Int_t i = 0; i < 3; i++)
361 target.fCutClusterRequirementITS[i] = fCutClusterRequirementITS[i];
363 target.fCutMaxC11 = fCutMaxC11;
364 target.fCutMaxC22 = fCutMaxC22;
365 target.fCutMaxC33 = fCutMaxC33;
366 target.fCutMaxC44 = fCutMaxC44;
367 target.fCutMaxC55 = fCutMaxC55;
369 target.fCutAcceptKinkDaughters = fCutAcceptKinkDaughters;
370 target.fCutRequireTPCRefit = fCutRequireTPCRefit;
371 target.fCutRequireITSRefit = fCutRequireITSRefit;
373 target.fCutNsigmaToVertex = fCutNsigmaToVertex;
374 target.fCutSigmaToVertexRequired = fCutSigmaToVertexRequired;
375 target.fCutMaxDCAToVertexXY = fCutMaxDCAToVertexXY;
376 target.fCutMaxDCAToVertexZ = fCutMaxDCAToVertexZ;
377 target.fCutDCAToVertex2D = fCutDCAToVertex2D;
378 target.fCutMinDCAToVertexXY = fCutMinDCAToVertexXY;
379 target.fCutMinDCAToVertexZ = fCutMinDCAToVertexZ;
381 target.fPMin = fPMin;
382 target.fPMax = fPMax;
383 target.fPtMin = fPtMin;
384 target.fPtMax = fPtMax;
385 target.fPxMin = fPxMin;
386 target.fPxMax = fPxMax;
387 target.fPyMin = fPyMin;
388 target.fPyMax = fPyMax;
389 target.fPzMin = fPzMin;
390 target.fPzMax = fPzMax;
391 target.fEtaMin = fEtaMin;
392 target.fEtaMax = fEtaMax;
393 target.fRapMin = fRapMin;
394 target.fRapMax = fRapMax;
396 target.fHistogramsOn = fHistogramsOn;
398 for (Int_t i=0; i<2; ++i)
400 if (fhNClustersITS[i]) target.fhNClustersITS[i] = (TH1F*) fhNClustersITS[i]->Clone();
401 if (fhNClustersTPC[i]) target.fhNClustersTPC[i] = (TH1F*) fhNClustersTPC[i]->Clone();
403 if (fhChi2PerClusterITS[i]) target.fhChi2PerClusterITS[i] = (TH1F*) fhChi2PerClusterITS[i]->Clone();
404 if (fhChi2PerClusterTPC[i]) target.fhChi2PerClusterTPC[i] = (TH1F*) fhChi2PerClusterTPC[i]->Clone();
406 if (fhC11[i]) target.fhC11[i] = (TH1F*) fhC11[i]->Clone();
407 if (fhC22[i]) target.fhC22[i] = (TH1F*) fhC22[i]->Clone();
408 if (fhC33[i]) target.fhC33[i] = (TH1F*) fhC33[i]->Clone();
409 if (fhC44[i]) target.fhC44[i] = (TH1F*) fhC44[i]->Clone();
410 if (fhC55[i]) target.fhC55[i] = (TH1F*) fhC55[i]->Clone();
412 if (fhDXY[i]) target.fhDXY[i] = (TH1F*) fhDXY[i]->Clone();
413 if (fhDZ[i]) target.fhDZ[i] = (TH1F*) fhDZ[i]->Clone();
414 if (fhDXYDZ[i]) target.fhDXYDZ[i] = (TH1F*) fhDXYDZ[i]->Clone();
415 if (fhDXYvsDZ[i]) target.fhDXYvsDZ[i] = (TH2F*) fhDXYvsDZ[i]->Clone();
417 if (fhDXYNormalized[i]) target.fhDXYNormalized[i] = (TH1F*) fhDXYNormalized[i]->Clone();
418 if (fhDZNormalized[i]) target.fhDZNormalized[i] = (TH1F*) fhDZNormalized[i]->Clone();
419 if (fhDXYvsDZNormalized[i]) target.fhDXYvsDZNormalized[i] = (TH2F*) fhDXYvsDZNormalized[i]->Clone();
420 if (fhNSigmaToVertex[i]) target.fhNSigmaToVertex[i] = (TH1F*) fhNSigmaToVertex[i]->Clone();
422 if (fhPt[i]) target.fhPt[i] = (TH1F*) fhPt[i]->Clone();
423 if (fhEta[i]) target.fhEta[i] = (TH1F*) fhEta[i]->Clone();
425 if (ffDTheoretical) target.ffDTheoretical = (TF1*) ffDTheoretical->Clone();
427 if (fhCutStatistics) target.fhCutStatistics = (TH1F*) fhCutStatistics->Clone();
428 if (fhCutCorrelation) target.fhCutCorrelation = (TH2F*) fhCutCorrelation->Clone();
433 //_____________________________________________________________________________
434 Long64_t AliESDtrackCuts::Merge(TCollection* list) {
435 // Merge a list of AliESDtrackCuts objects with this (needed for PROOF)
436 // Returns the number of merged objects (including this)
443 TIterator* iter = list->MakeIterator();
446 // collection of measured and generated histograms
448 while ((obj = iter->Next())) {
450 AliESDtrackCuts* entry = dynamic_cast<AliESDtrackCuts*>(obj);
454 if (!entry->fHistogramsOn)
457 for (Int_t i=0; i<2; i++) {
459 fhNClustersITS[i] ->Add(entry->fhNClustersITS[i] );
460 fhNClustersTPC[i] ->Add(entry->fhNClustersTPC[i] );
462 fhChi2PerClusterITS[i] ->Add(entry->fhChi2PerClusterITS[i]);
463 fhChi2PerClusterTPC[i] ->Add(entry->fhChi2PerClusterTPC[i]);
465 fhC11[i] ->Add(entry->fhC11[i] );
466 fhC22[i] ->Add(entry->fhC22[i] );
467 fhC33[i] ->Add(entry->fhC33[i] );
468 fhC44[i] ->Add(entry->fhC44[i] );
469 fhC55[i] ->Add(entry->fhC55[i] );
471 fhDXY[i] ->Add(entry->fhDXY[i] );
472 fhDZ[i] ->Add(entry->fhDZ[i] );
473 fhDXYDZ[i] ->Add(entry->fhDXYDZ[i] );
474 fhDXYvsDZ[i] ->Add(entry->fhDXYvsDZ[i] );
476 fhDXYNormalized[i] ->Add(entry->fhDXYNormalized[i] );
477 fhDZNormalized[i] ->Add(entry->fhDZNormalized[i] );
478 fhDXYvsDZNormalized[i] ->Add(entry->fhDXYvsDZNormalized[i]);
479 fhNSigmaToVertex[i] ->Add(entry->fhNSigmaToVertex[i]);
481 fhPt[i] ->Add(entry->fhPt[i]);
482 fhEta[i] ->Add(entry->fhEta[i]);
485 fhCutStatistics ->Add(entry->fhCutStatistics);
486 fhCutCorrelation ->Add(entry->fhCutCorrelation);
494 //____________________________________________________________________
495 Float_t AliESDtrackCuts::GetSigmaToVertex(AliESDtrack* esdTrack)
497 // Calculates the number of sigma to the vertex.
502 esdTrack->GetImpactParameters(b,bCov);
504 if (bCov[0]<=0 || bCov[2]<=0) {
505 AliDebugClass(1, "Estimated b resolution lower or equal zero!");
506 bCov[0]=0; bCov[2]=0;
508 bRes[0] = TMath::Sqrt(bCov[0]);
509 bRes[1] = TMath::Sqrt(bCov[2]);
511 // -----------------------------------
512 // How to get to a n-sigma cut?
514 // The accumulated statistics from 0 to d is
516 // -> Erf(d/Sqrt(2)) for a 1-dim gauss (d = n_sigma)
517 // -> 1 - Exp(-d**2) for a 2-dim gauss (d*d = dx*dx + dy*dy != n_sigma)
519 // It means that for a 2-dim gauss: n_sigma(d) = Sqrt(2)*ErfInv(1 - Exp((-d**2)/2)
520 // Can this be expressed in a different way?
522 if (bRes[0] == 0 || bRes[1] ==0)
525 Float_t d = TMath::Sqrt(TMath::Power(b[0]/bRes[0],2) + TMath::Power(b[1]/bRes[1],2));
527 // work around precision problem
528 // if d is too big, TMath::Exp(...) gets 0, and TMath::ErfInverse(1) that should be infinite, gets 0 :(
529 // 1e-15 corresponds to nsigma ~ 7.7
530 if (TMath::Exp(-d * d / 2) < 1e-15)
533 Float_t nSigma = TMath::ErfInverse(1 - TMath::Exp(-d * d / 2)) * TMath::Sqrt(2);
537 void AliESDtrackCuts::EnableNeededBranches(TTree* tree)
539 // enables the branches needed by AcceptTrack, for a list see comment of AcceptTrack
541 tree->SetBranchStatus("fTracks.fFlags", 1);
542 tree->SetBranchStatus("fTracks.fITSncls", 1);
543 tree->SetBranchStatus("fTracks.fTPCncls", 1);
544 tree->SetBranchStatus("fTracks.fITSchi2", 1);
545 tree->SetBranchStatus("fTracks.fTPCchi2", 1);
546 tree->SetBranchStatus("fTracks.fC*", 1);
547 tree->SetBranchStatus("fTracks.fD", 1);
548 tree->SetBranchStatus("fTracks.fZ", 1);
549 tree->SetBranchStatus("fTracks.fCdd", 1);
550 tree->SetBranchStatus("fTracks.fCdz", 1);
551 tree->SetBranchStatus("fTracks.fCzz", 1);
552 tree->SetBranchStatus("fTracks.fP*", 1);
553 tree->SetBranchStatus("fTracks.fR*", 1);
554 tree->SetBranchStatus("fTracks.fKinkIndexes*", 1);
557 //____________________________________________________________________
559 AliESDtrackCuts::AcceptTrack(AliESDtrack* esdTrack) {
561 // figure out if the tracks survives all the track cuts defined
563 // the different quality parameter and kinematic values are first
564 // retrieved from the track. then it is found out what cuts the
565 // track did not survive and finally the cuts are imposed.
567 // this function needs the following branches:
573 // fTracks.fC //GetExternalCovariance
574 // fTracks.fD //GetImpactParameters
575 // fTracks.fZ //GetImpactParameters
576 // fTracks.fCdd //GetImpactParameters
577 // fTracks.fCdz //GetImpactParameters
578 // fTracks.fCzz //GetImpactParameters
579 // fTracks.fP //GetPxPyPz
580 // fTracks.fR //GetMass
581 // fTracks.fP //GetMass
582 // fTracks.fKinkIndexes
584 UInt_t status = esdTrack->GetStatus();
586 // getting quality parameters from the ESD track
587 Int_t nClustersITS = esdTrack->GetITSclusters(0);
588 Int_t nClustersTPC = esdTrack->GetTPCclusters(0);
590 Float_t chi2PerClusterITS = -1;
591 Float_t chi2PerClusterTPC = -1;
593 chi2PerClusterITS = esdTrack->GetITSchi2()/Float_t(nClustersITS);
595 chi2PerClusterTPC = esdTrack->GetTPCchi2()/Float_t(nClustersTPC);
597 esdTrack->GetExternalCovariance(extCov);
599 // getting the track to vertex parameters
600 Float_t nSigmaToVertex = GetSigmaToVertex(esdTrack);
604 esdTrack->GetImpactParameters(b,bCov);
605 if (bCov[0]<=0 || bCov[2]<=0) {
606 AliDebug(1, "Estimated b resolution lower or equal zero!");
607 bCov[0]=0; bCov[2]=0;
610 Float_t dcaToVertexXY = b[0];
611 Float_t dcaToVertexZ = b[1];
613 Float_t dcaToVertex = -1;
615 if (fCutDCAToVertex2D)
617 dcaToVertex = TMath::Sqrt(dcaToVertexXY*dcaToVertexXY/fCutMaxDCAToVertexXY/fCutMaxDCAToVertexXY + dcaToVertexZ*dcaToVertexZ/fCutMaxDCAToVertexZ/fCutMaxDCAToVertexZ);
620 dcaToVertex = TMath::Sqrt(dcaToVertexXY*dcaToVertexXY + dcaToVertexZ*dcaToVertexZ);
622 // getting the kinematic variables of the track
623 // (assuming the mass is known)
625 esdTrack->GetPxPyPz(p);
627 Float_t momentum = TMath::Sqrt(TMath::Power(p[0],2) + TMath::Power(p[1],2) + TMath::Power(p[2],2));
628 Float_t pt = TMath::Sqrt(TMath::Power(p[0],2) + TMath::Power(p[1],2));
629 Float_t energy = TMath::Sqrt(TMath::Power(esdTrack->GetMass(),2) + TMath::Power(momentum,2));
632 //y-eta related calculations
635 if((momentum != TMath::Abs(p[2]))&&(momentum != 0))
636 eta = 0.5*TMath::Log((momentum + p[2])/(momentum - p[2]));
637 if((energy != TMath::Abs(p[2]))&&(momentum != 0))
638 y = 0.5*TMath::Log((energy + p[2])/(energy - p[2]));
641 //########################################################################
645 for (Int_t i=0; i<kNCuts; i++) cuts[i]=kFALSE;
647 // track quality cuts
648 if (fCutRequireTPCRefit && (status&AliESDtrack::kTPCrefit)==0)
650 if (fCutRequireITSRefit && (status&AliESDtrack::kITSrefit)==0)
652 if (nClustersTPC<fCutMinNClusterTPC)
654 if (nClustersITS<fCutMinNClusterITS)
656 if (chi2PerClusterTPC>fCutMaxChi2PerClusterTPC)
658 if (chi2PerClusterITS>fCutMaxChi2PerClusterITS)
660 if (extCov[0] > fCutMaxC11)
662 if (extCov[2] > fCutMaxC22)
664 if (extCov[5] > fCutMaxC33)
666 if (extCov[9] > fCutMaxC44)
668 if (extCov[14] > fCutMaxC55)
670 if (nSigmaToVertex > fCutNsigmaToVertex && fCutSigmaToVertexRequired)
672 // if n sigma could not be calculated
673 if (nSigmaToVertex<0 && fCutSigmaToVertexRequired)
675 if (!fCutAcceptKinkDaughters && esdTrack->GetKinkIndex(0)>0)
677 // track kinematics cut
678 if((momentum < fPMin) || (momentum > fPMax))
680 if((pt < fPtMin) || (pt > fPtMax))
682 if((p[0] < fPxMin) || (p[0] > fPxMax))
684 if((p[1] < fPyMin) || (p[1] > fPyMax))
686 if((p[2] < fPzMin) || (p[2] > fPzMax))
688 if((eta < fEtaMin) || (eta > fEtaMax))
690 if((y < fRapMin) || (y > fRapMax))
692 if (fCutDCAToVertex2D && dcaToVertex > 1)
694 if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexXY) > fCutMaxDCAToVertexXY)
696 if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexZ) > fCutMaxDCAToVertexZ)
698 if (fCutDCAToVertex2D && fCutMinDCAToVertexXY > 0 && fCutMinDCAToVertexZ > 0 && dcaToVertexXY*dcaToVertexXY/fCutMinDCAToVertexXY/fCutMinDCAToVertexXY + dcaToVertexZ*dcaToVertexZ/fCutMinDCAToVertexZ/fCutMinDCAToVertexZ < 1)
700 if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexXY) < fCutMinDCAToVertexXY)
702 if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexZ) < fCutMinDCAToVertexZ)
705 for (Int_t i = 0; i < 3; i++)
706 cuts[27+i] = !CheckITSClusterRequirement(fCutClusterRequirementITS[i], esdTrack->HasPointOnITSLayer(i*2), esdTrack->HasPointOnITSLayer(i*2+1));
709 for (Int_t i=0; i<kNCuts; i++)
710 if (cuts[i]) {cut = kTRUE;}
714 //########################################################################
715 // filling histograms
717 fhCutStatistics->Fill(fhCutStatistics->GetBinCenter(fhCutStatistics->GetXaxis()->FindBin("n tracks")));
719 fhCutStatistics->Fill(fhCutStatistics->GetBinCenter(fhCutStatistics->GetXaxis()->FindBin("n cut tracks")));
721 for (Int_t i=0; i<kNCuts; i++) {
723 fhCutStatistics->Fill(fhCutStatistics->GetBinCenter(fhCutStatistics->GetXaxis()->FindBin(fgkCutNames[i])));
725 for (Int_t j=i; j<kNCuts; j++) {
726 if (cuts[i] && cuts[j]) {
727 Float_t xC = fhCutCorrelation->GetXaxis()->GetBinCenter(fhCutCorrelation->GetXaxis()->FindBin(fgkCutNames[i]));
728 Float_t yC = fhCutCorrelation->GetYaxis()->GetBinCenter(fhCutCorrelation->GetYaxis()->FindBin(fgkCutNames[j]));
729 fhCutCorrelation->Fill(xC, yC);
735 // now we loop over the filling of the histograms twice: once "before" the cut, once "after"
736 // the code is not in a function due to too many local variables that would need to be passed
738 for (Int_t id = 0; id < 2; id++)
740 // id = 0 --> before cut
741 // id = 1 --> after cut
745 fhNClustersITS[id]->Fill(nClustersITS);
746 fhNClustersTPC[id]->Fill(nClustersTPC);
747 fhChi2PerClusterITS[id]->Fill(chi2PerClusterITS);
748 fhChi2PerClusterTPC[id]->Fill(chi2PerClusterTPC);
750 fhC11[id]->Fill(extCov[0]);
751 fhC22[id]->Fill(extCov[2]);
752 fhC33[id]->Fill(extCov[5]);
753 fhC44[id]->Fill(extCov[9]);
754 fhC55[id]->Fill(extCov[14]);
757 fhEta[id]->Fill(eta);
760 bRes[0] = TMath::Sqrt(bCov[0]);
761 bRes[1] = TMath::Sqrt(bCov[2]);
763 fhDZ[id]->Fill(b[1]);
764 fhDXY[id]->Fill(b[0]);
765 fhDXYDZ[id]->Fill(dcaToVertex);
766 fhDXYvsDZ[id]->Fill(b[1],b[0]);
768 if (bRes[0]!=0 && bRes[1]!=0) {
769 fhDZNormalized[id]->Fill(b[1]/bRes[1]);
770 fhDXYNormalized[id]->Fill(b[0]/bRes[0]);
771 fhDXYvsDZNormalized[id]->Fill(b[1]/bRes[1], b[0]/bRes[0]);
772 fhNSigmaToVertex[id]->Fill(nSigmaToVertex);
784 //____________________________________________________________________
785 Bool_t AliESDtrackCuts::CheckITSClusterRequirement(ITSClusterRequirement req, Bool_t clusterL1, Bool_t clusterL2)
787 // checks if the cluster requirement is fullfilled (in this case: return kTRUE)
791 case kOff: return kTRUE;
792 case kNone: return !clusterL1 && !clusterL2;
793 case kAny: return clusterL1 || clusterL2;
794 case kFirst: return clusterL1;
795 case kOnlyFirst: return clusterL1 && !clusterL2;
796 case kSecond: return clusterL2;
797 case kOnlySecond: return clusterL2 && !clusterL1;
798 case kBoth: return clusterL1 && clusterL2;
804 //____________________________________________________________________
805 AliESDtrack* AliESDtrackCuts::GetTPCOnlyTrack(AliESDEvent* esd, Int_t iTrack)
807 // creates a TPC only track from the given esd track
808 // the track has to be deleted by the user
810 // NB. most of the functionality to get a TPC only track from an ESD track is in AliESDtrack, where it should be
811 // there are only missing propagations here that are needed for old data
812 // this function will therefore become obsolete
814 // adapted from code provided by CKB
816 if (!esd->GetPrimaryVertexTPC())
817 return 0; // No TPC vertex no TPC tracks
819 if(!esd->GetPrimaryVertexTPC()->GetStatus())
820 return 0; // TPC Vertex is created by default in AliESDEvent, do not use in this case
822 AliESDtrack* track = esd->GetTrack(iTrack);
826 AliESDtrack *tpcTrack = new AliESDtrack();
828 // This should have been done during the reconstruction
829 // fixed by Juri in r26675
830 // but recalculate for older data CKB
832 track->GetImpactParametersTPC(p,cov);
834 track->RelateToVertexTPC(esd->GetPrimaryVertexTPC(),esd->GetMagneticField(),kVeryBig);
837 // only true if we have a tpc track
838 if (!track->FillTPCOnlyTrack(*tpcTrack))
844 // propagate to Vertex
845 // not needed for normal reconstructed ESDs...
846 // Double_t pTPC[2],covTPC[3];
847 // tpcTrack->PropagateToDCA(esd->GetPrimaryVertexTPC(), esd->GetMagneticField(), 10000, pTPC, covTPC);
852 //____________________________________________________________________
853 TObjArray* AliESDtrackCuts::GetAcceptedTracks(AliESDEvent* esd,Bool_t bTPC)
856 // returns an array of all tracks that pass the cuts
857 // or an array of TPC only tracks (propagated to the TPC vertex during reco)
858 // tracks that pass the cut
860 TObjArray* acceptedTracks = new TObjArray();
862 // loop over esd tracks
863 for (Int_t iTrack = 0; iTrack < esd->GetNumberOfTracks(); iTrack++) {
865 if(!esd->GetPrimaryVertexTPC())return acceptedTracks; // No TPC vertex no TPC tracks
866 if(!esd->GetPrimaryVertexTPC()->GetStatus())return acceptedTracks; // No proper TPC vertex, only the default
868 AliESDtrack *tpcTrack = GetTPCOnlyTrack(esd, iTrack);
872 if (AcceptTrack(tpcTrack)) {
873 acceptedTracks->Add(tpcTrack);
880 AliESDtrack* track = esd->GetTrack(iTrack);
881 if(AcceptTrack(track))
882 acceptedTracks->Add(track);
885 if(bTPC)acceptedTracks->SetOwner(kTRUE);
886 return acceptedTracks;
889 //____________________________________________________________________
890 Int_t AliESDtrackCuts::CountAcceptedTracks(AliESDEvent* esd)
893 // returns an the number of tracks that pass the cuts
898 // loop over esd tracks
899 for (Int_t iTrack = 0; iTrack < esd->GetNumberOfTracks(); iTrack++) {
900 AliESDtrack* track = esd->GetTrack(iTrack);
901 if (AcceptTrack(track))
908 //____________________________________________________________________
909 void AliESDtrackCuts::DefineHistograms(Int_t color) {
911 // diagnostics histograms are defined
916 Bool_t oldStatus = TH1::AddDirectoryStatus();
917 TH1::AddDirectory(kFALSE);
919 //###################################################################################
920 // defining histograms
922 fhCutStatistics = new TH1F("cut_statistics","cut statistics",kNCuts+4,-0.5,kNCuts+3.5);
924 fhCutStatistics->GetXaxis()->SetBinLabel(1,"n tracks");
925 fhCutStatistics->GetXaxis()->SetBinLabel(2,"n cut tracks");
927 fhCutCorrelation = new TH2F("cut_correlation","cut correlation",kNCuts,-0.5,kNCuts-0.5,kNCuts,-0.5,kNCuts-0.5);;
929 for (Int_t i=0; i<kNCuts; i++) {
930 fhCutStatistics->GetXaxis()->SetBinLabel(i+4,fgkCutNames[i]);
931 fhCutCorrelation->GetXaxis()->SetBinLabel(i+1,fgkCutNames[i]);
932 fhCutCorrelation->GetYaxis()->SetBinLabel(i+1,fgkCutNames[i]);
935 fhCutStatistics ->SetLineColor(color);
936 fhCutCorrelation ->SetLineColor(color);
937 fhCutStatistics ->SetLineWidth(2);
938 fhCutCorrelation ->SetLineWidth(2);
940 for (Int_t i=0; i<2; i++) {
941 fhNClustersITS[i] = new TH1F("nClustersITS" ,"",8,-0.5,7.5);
942 fhNClustersTPC[i] = new TH1F("nClustersTPC" ,"",165,-0.5,164.5);
943 fhChi2PerClusterITS[i] = new TH1F("chi2PerClusterITS","",500,0,10);
944 fhChi2PerClusterTPC[i] = new TH1F("chi2PerClusterTPC","",500,0,10);
946 fhC11[i] = new TH1F("covMatrixDiagonal11","",2000,0,20);
947 fhC22[i] = new TH1F("covMatrixDiagonal22","",2000,0,20);
948 fhC33[i] = new TH1F("covMatrixDiagonal33","",1000,0,0.1);
949 fhC44[i] = new TH1F("covMatrixDiagonal44","",1000,0,0.1);
950 fhC55[i] = new TH1F("covMatrixDiagonal55","",1000,0,5);
952 fhDXY[i] = new TH1F("dXY" ,"",500,-10,10);
953 fhDZ[i] = new TH1F("dZ" ,"",500,-10,10);
954 fhDXYDZ[i] = new TH1F("dXYDZ" ,"",500,0,10);
955 fhDXYvsDZ[i] = new TH2F("dXYvsDZ","",200,-10,10,200,-10,10);
957 fhDXYNormalized[i] = new TH1F("dXYNormalized" ,"",500,-10,10);
958 fhDZNormalized[i] = new TH1F("dZNormalized" ,"",500,-10,10);
959 fhDXYvsDZNormalized[i] = new TH2F("dXYvsDZNormalized","",200,-10,10,200,-10,10);
961 fhNSigmaToVertex[i] = new TH1F("nSigmaToVertex","",500,0,10);
963 fhPt[i] = new TH1F("pt" ,"p_{T} distribution;p_{T} (GeV/c)", 800, 0.0, 10.0);
964 fhEta[i] = new TH1F("eta" ,"#eta distribution;#eta",40,-2.0,2.0);
966 fhNClustersITS[i]->SetTitle("n ITS clusters");
967 fhNClustersTPC[i]->SetTitle("n TPC clusters");
968 fhChi2PerClusterITS[i]->SetTitle("#Chi^{2} per ITS cluster");
969 fhChi2PerClusterTPC[i]->SetTitle("#Chi^{2} per TPC cluster");
971 fhC11[i]->SetTitle("cov 11 : #sigma_{y}^{2} [cm^{2}]");
972 fhC22[i]->SetTitle("cov 22 : #sigma_{z}^{2} [cm^{2}]");
973 fhC33[i]->SetTitle("cov 33 : #sigma_{sin(#phi)}^{2}");
974 fhC44[i]->SetTitle("cov 44 : #sigma_{tan(#theta_{dip})}^{2}");
975 fhC55[i]->SetTitle("cov 55 : #sigma_{1/p_{T}}^{2} [(c/GeV)^2]");
977 fhDXY[i]->SetXTitle("transverse impact parameter (cm)");
978 fhDZ[i]->SetXTitle("longitudinal impact parameter (cm)");
979 fhDXYDZ[i]->SetTitle("absolute impact parameter;sqrt(dXY**2 + dZ**2) (cm)");
980 fhDXYvsDZ[i]->SetXTitle("longitudinal impact parameter (cm)");
981 fhDXYvsDZ[i]->SetYTitle("transverse impact parameter (cm)");
983 fhDXYNormalized[i]->SetTitle("normalized trans impact par (n#sigma)");
984 fhDZNormalized[i]->SetTitle("normalized long impact par (n#sigma)");
985 fhDXYvsDZNormalized[i]->SetTitle("normalized long impact par (n#sigma)");
986 fhDXYvsDZNormalized[i]->SetYTitle("normalized trans impact par (n#sigma)");
987 fhNSigmaToVertex[i]->SetTitle("n #sigma to vertex");
989 fhNClustersITS[i]->SetLineColor(color); fhNClustersITS[i]->SetLineWidth(2);
990 fhNClustersTPC[i]->SetLineColor(color); fhNClustersTPC[i]->SetLineWidth(2);
991 fhChi2PerClusterITS[i]->SetLineColor(color); fhChi2PerClusterITS[i]->SetLineWidth(2);
992 fhChi2PerClusterTPC[i]->SetLineColor(color); fhChi2PerClusterTPC[i]->SetLineWidth(2);
994 fhC11[i]->SetLineColor(color); fhC11[i]->SetLineWidth(2);
995 fhC22[i]->SetLineColor(color); fhC22[i]->SetLineWidth(2);
996 fhC33[i]->SetLineColor(color); fhC33[i]->SetLineWidth(2);
997 fhC44[i]->SetLineColor(color); fhC44[i]->SetLineWidth(2);
998 fhC55[i]->SetLineColor(color); fhC55[i]->SetLineWidth(2);
1000 fhDXY[i]->SetLineColor(color); fhDXY[i]->SetLineWidth(2);
1001 fhDZ[i]->SetLineColor(color); fhDZ[i]->SetLineWidth(2);
1002 fhDXYDZ[i]->SetLineColor(color); fhDXYDZ[i]->SetLineWidth(2);
1004 fhDXYNormalized[i]->SetLineColor(color); fhDXYNormalized[i]->SetLineWidth(2);
1005 fhDZNormalized[i]->SetLineColor(color); fhDZNormalized[i]->SetLineWidth(2);
1006 fhNSigmaToVertex[i]->SetLineColor(color); fhNSigmaToVertex[i]->SetLineWidth(2);
1009 // The number of sigmas to the vertex is per definition gaussian
1010 ffDTheoretical = new TF1("nSigmaToVertexTheoretical","([0]/2.506628274)*exp(-(x**2)/2)",0,50);
1011 ffDTheoretical->SetParameter(0,1);
1013 TH1::AddDirectory(oldStatus);
1016 //____________________________________________________________________
1017 Bool_t AliESDtrackCuts::LoadHistograms(const Char_t* dir)
1020 // loads the histograms from a file
1021 // if dir is empty a directory with the name of this object is taken (like in SaveHistogram)
1027 if (!gDirectory->cd(dir))
1030 ffDTheoretical = dynamic_cast<TF1*> (gDirectory->Get("nSigmaToVertexTheory"));
1032 fhCutStatistics = dynamic_cast<TH1F*> (gDirectory->Get("cut_statistics"));
1033 fhCutCorrelation = dynamic_cast<TH2F*> (gDirectory->Get("cut_correlation"));
1035 for (Int_t i=0; i<2; i++) {
1038 gDirectory->cd("before_cuts");
1041 gDirectory->cd("after_cuts");
1043 fhNClustersITS[i] = dynamic_cast<TH1F*> (gDirectory->Get("nClustersITS" ));
1044 fhNClustersTPC[i] = dynamic_cast<TH1F*> (gDirectory->Get("nClustersTPC" ));
1045 fhChi2PerClusterITS[i] = dynamic_cast<TH1F*> (gDirectory->Get("chi2PerClusterITS"));
1046 fhChi2PerClusterTPC[i] = dynamic_cast<TH1F*> (gDirectory->Get("chi2PerClusterTPC"));
1048 fhC11[i] = dynamic_cast<TH1F*> (gDirectory->Get("covMatrixDiagonal11"));
1049 fhC22[i] = dynamic_cast<TH1F*> (gDirectory->Get("covMatrixDiagonal22"));
1050 fhC33[i] = dynamic_cast<TH1F*> (gDirectory->Get("covMatrixDiagonal33"));
1051 fhC44[i] = dynamic_cast<TH1F*> (gDirectory->Get("covMatrixDiagonal44"));
1052 fhC55[i] = dynamic_cast<TH1F*> (gDirectory->Get("covMatrixDiagonal55"));
1054 fhDXY[i] = dynamic_cast<TH1F*> (gDirectory->Get("dXY" ));
1055 fhDZ[i] = dynamic_cast<TH1F*> (gDirectory->Get("dZ" ));
1056 fhDXYDZ[i] = dynamic_cast<TH1F*> (gDirectory->Get("dXYDZ"));
1057 fhDXYvsDZ[i] = dynamic_cast<TH2F*> (gDirectory->Get("dXYvsDZ"));
1059 fhDXYNormalized[i] = dynamic_cast<TH1F*> (gDirectory->Get("dXYNormalized" ));
1060 fhDZNormalized[i] = dynamic_cast<TH1F*> (gDirectory->Get("dZNormalized" ));
1061 fhDXYvsDZNormalized[i] = dynamic_cast<TH2F*> (gDirectory->Get("dXYvsDZNormalized"));
1062 fhNSigmaToVertex[i] = dynamic_cast<TH1F*> (gDirectory->Get("nSigmaToVertex"));
1064 fhPt[i] = dynamic_cast<TH1F*> (gDirectory->Get("pt"));
1065 fhEta[i] = dynamic_cast<TH1F*> (gDirectory->Get("eta"));
1067 gDirectory->cd("../");
1070 gDirectory->cd("..");
1075 //____________________________________________________________________
1076 void AliESDtrackCuts::SaveHistograms(const Char_t* dir) {
1078 // saves the histograms in a directory (dir)
1081 if (!fHistogramsOn) {
1082 AliDebug(0, "Histograms not on - cannot save histograms!!!");
1089 gDirectory->mkdir(dir);
1090 gDirectory->cd(dir);
1092 gDirectory->mkdir("before_cuts");
1093 gDirectory->mkdir("after_cuts");
1095 // a factor of 2 is needed since n sigma is positive
1096 ffDTheoretical->SetParameter(0,2*fhNSigmaToVertex[0]->Integral("width"));
1097 ffDTheoretical->Write("nSigmaToVertexTheory");
1099 fhCutStatistics->Write();
1100 fhCutCorrelation->Write();
1102 for (Int_t i=0; i<2; i++) {
1104 gDirectory->cd("before_cuts");
1106 gDirectory->cd("after_cuts");
1108 fhNClustersITS[i] ->Write();
1109 fhNClustersTPC[i] ->Write();
1110 fhChi2PerClusterITS[i] ->Write();
1111 fhChi2PerClusterTPC[i] ->Write();
1121 fhDXYDZ[i] ->Write();
1122 fhDXYvsDZ[i] ->Write();
1124 fhDXYNormalized[i] ->Write();
1125 fhDZNormalized[i] ->Write();
1126 fhDXYvsDZNormalized[i] ->Write();
1127 fhNSigmaToVertex[i] ->Write();
1132 gDirectory->cd("../");
1135 gDirectory->cd("../");
1138 //____________________________________________________________________
1139 void AliESDtrackCuts::DrawHistograms()
1141 // draws some histograms
1143 TCanvas* canvas1 = new TCanvas(Form("%s_1", GetName()), "Track Quality Results1", 800, 800);
1144 canvas1->Divide(2, 2);
1147 fhNClustersTPC[0]->SetStats(kFALSE);
1148 fhNClustersTPC[0]->Draw();
1151 fhChi2PerClusterTPC[0]->SetStats(kFALSE);
1152 fhChi2PerClusterTPC[0]->Draw();
1155 fhNSigmaToVertex[0]->SetStats(kFALSE);
1156 fhNSigmaToVertex[0]->GetXaxis()->SetRangeUser(0, 10);
1157 fhNSigmaToVertex[0]->Draw();
1159 canvas1->SaveAs(Form("%s_%s.gif", GetName(), canvas1->GetName()));
1161 TCanvas* canvas2 = new TCanvas(Form("%s_2", GetName()), "Track Quality Results2", 1200, 800);
1162 canvas2->Divide(3, 2);
1165 fhC11[0]->SetStats(kFALSE);
1170 fhC22[0]->SetStats(kFALSE);
1175 fhC33[0]->SetStats(kFALSE);
1180 fhC44[0]->SetStats(kFALSE);
1185 fhC55[0]->SetStats(kFALSE);
1189 canvas2->SaveAs(Form("%s_%s.gif", GetName(), canvas2->GetName()));
1191 TCanvas* canvas3 = new TCanvas(Form("%s_3", GetName()), "Track Quality Results3", 1200, 800);
1192 canvas3->Divide(3, 2);
1195 fhDXY[0]->SetStats(kFALSE);
1200 fhDZ[0]->SetStats(kFALSE);
1205 fhDXYvsDZ[0]->SetStats(kFALSE);
1207 gPad->SetRightMargin(0.15);
1208 fhDXYvsDZ[0]->Draw("COLZ");
1211 fhDXYNormalized[0]->SetStats(kFALSE);
1213 fhDXYNormalized[0]->Draw();
1216 fhDZNormalized[0]->SetStats(kFALSE);
1218 fhDZNormalized[0]->Draw();
1221 fhDXYvsDZNormalized[0]->SetStats(kFALSE);
1223 gPad->SetRightMargin(0.15);
1224 fhDXYvsDZNormalized[0]->Draw("COLZ");
1226 canvas3->SaveAs(Form("%s_%s.gif", GetName(), canvas3->GetName()));
1228 TCanvas* canvas4 = new TCanvas(Form("%s_4", GetName()), "Track Quality Results4", 800, 500);
1229 canvas4->Divide(2, 1);
1232 fhCutStatistics->SetStats(kFALSE);
1233 fhCutStatistics->LabelsOption("v");
1234 gPad->SetBottomMargin(0.3);
1235 fhCutStatistics->Draw();
1238 fhCutCorrelation->SetStats(kFALSE);
1239 fhCutCorrelation->LabelsOption("v");
1240 gPad->SetBottomMargin(0.3);
1241 gPad->SetLeftMargin(0.3);
1242 fhCutCorrelation->Draw("COLZ");
1244 canvas4->SaveAs(Form("%s_%s.gif", GetName(), canvas4->GetName()));
1247 fhDXYvsDZNormalized[0]->SetStats(kFALSE);
1248 fhDXYvsDZNormalized[0]->DrawCopy("COLZ");
1251 fhNClustersTPC[0]->SetStats(kFALSE);
1252 fhNClustersTPC[0]->DrawCopy();
1255 fhChi2PerClusterITS[0]->SetStats(kFALSE);
1256 fhChi2PerClusterITS[0]->DrawCopy();
1257 fhChi2PerClusterITS[1]->SetLineColor(2);
1258 fhChi2PerClusterITS[1]->DrawCopy("SAME");
1261 fhChi2PerClusterTPC[0]->SetStats(kFALSE);
1262 fhChi2PerClusterTPC[0]->DrawCopy();
1263 fhChi2PerClusterTPC[1]->SetLineColor(2);
1264 fhChi2PerClusterTPC[1]->DrawCopy("SAME");*/
1267 Float_t AliESDtrackCuts::GetMinNsigmaToVertex() const
1269 // deprecated, please use GetMaxNsigmaToVertex
1271 Printf("WARNING: AliESDtrackCuts::GetMinNsigmaToVertex is DEPRECATED and will be removed in the next release. Please use GetMaxNsigmaToVertex instead. Renaming was done to improve code readability.");
1273 return GetMaxNsigmaToVertex();
1276 void AliESDtrackCuts::SetMinNsigmaToVertex(Float_t sigma)
1278 // deprecated, will be removed in next release
1280 Printf("WARNING: AliESDtrackCuts::SetMinNsigmaToVertex is DEPRECATED and will be removed in the next release. Please use SetMaxNsigmaToVertex instead. Renaming was done to improve code readability.");
1282 SetMaxNsigmaToVertex(sigma);
1285 void AliESDtrackCuts::SetAcceptKingDaughters(Bool_t b)
1287 // deprecated, will be removed in next release
1289 Printf("WARNING: AliESDtrackCuts::SetAcceptKingDaughters is DEPRECATED and will be removed in the next release. Please use SetAcceptKinkDaughters instead. Renaming was done to improve code readability.");
1291 SetAcceptKinkDaughters(b);
1294 Bool_t AliESDtrackCuts::GetAcceptKingDaughters() const
1296 // deprecated, will be removed in next release
1298 Printf("WARNING: AliESDtrackCuts::GetAcceptKingDaughters is DEPRECATED and will be removed in the next release. Please use GetAcceptKinkDaughters instead. Renaming was done to improve code readability.");
1300 return GetAcceptKinkDaughters();