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 **************************************************************************/
18 #include "AliESDtrackCuts.h"
20 #include <AliESDtrack.h>
22 #include <AliESDEvent.h>
27 #include <TDirectory.h>
29 //____________________________________________________________________
30 ClassImp(AliESDtrackCuts)
33 const Char_t* AliESDtrackCuts::fgkCutNames[kNCuts] = {
38 "#Chi^{2}/clusters TPC",
39 "#Chi^{2}/clusters ITS",
57 //____________________________________________________________________
58 AliESDtrackCuts::AliESDtrackCuts(const Char_t* name, const Char_t* title) : AliAnalysisCuts(name,title),
59 fCutMinNClusterTPC(0),
60 fCutMinNClusterITS(0),
61 fCutMaxChi2PerClusterTPC(0),
62 fCutMaxChi2PerClusterITS(0),
68 fCutAcceptKinkDaughters(0),
69 fCutRequireTPCRefit(0),
70 fCutRequireITSRefit(0),
71 fCutNsigmaToVertex(0),
72 fCutSigmaToVertexRequired(0),
97 //##############################################################################
98 // setting default cuts
100 SetMinNClustersITS();
101 SetMaxChi2PerClusterTPC();
102 SetMaxChi2PerClusterITS();
103 SetMaxCovDiagonalElements();
104 SetRequireTPCRefit();
105 SetRequireITSRefit();
106 SetAcceptKingDaughters();
107 SetMinNsigmaToVertex();
108 SetRequireSigmaToVertex();
120 //_____________________________________________________________________________
121 AliESDtrackCuts::AliESDtrackCuts(const AliESDtrackCuts &c) : AliAnalysisCuts(c),
122 fCutMinNClusterTPC(0),
123 fCutMinNClusterITS(0),
124 fCutMaxChi2PerClusterTPC(0),
125 fCutMaxChi2PerClusterITS(0),
131 fCutAcceptKinkDaughters(0),
132 fCutRequireTPCRefit(0),
133 fCutRequireITSRefit(0),
134 fCutNsigmaToVertex(0),
135 fCutSigmaToVertexRequired(0),
159 ((AliESDtrackCuts &) c).Copy(*this);
162 AliESDtrackCuts::~AliESDtrackCuts()
168 for (Int_t i=0; i<2; i++) {
170 if (fhNClustersITS[i])
171 delete fhNClustersITS[i];
172 if (fhNClustersTPC[i])
173 delete fhNClustersTPC[i];
174 if (fhChi2PerClusterITS[i])
175 delete fhChi2PerClusterITS[i];
176 if (fhChi2PerClusterTPC[i])
177 delete fhChi2PerClusterTPC[i];
196 if (fhDXYNormalized[i])
197 delete fhDXYNormalized[i];
198 if (fhDZNormalized[i])
199 delete fhDZNormalized[i];
200 if (fhDXYvsDZNormalized[i])
201 delete fhDXYvsDZNormalized[i];
202 if (fhNSigmaToVertex[i])
203 delete fhNSigmaToVertex[i];
207 delete ffDTheoretical;
210 delete fhCutStatistics;
211 if (fhCutCorrelation)
212 delete fhCutCorrelation;
215 void AliESDtrackCuts::Init()
218 // sets everything to zero
221 fCutMinNClusterTPC = 0;
222 fCutMinNClusterITS = 0;
224 fCutMaxChi2PerClusterTPC = 0;
225 fCutMaxChi2PerClusterITS = 0;
233 fCutAcceptKinkDaughters = 0;
234 fCutRequireTPCRefit = 0;
235 fCutRequireITSRefit = 0;
237 fCutNsigmaToVertex = 0;
238 fCutSigmaToVertexRequired = 0;
255 fHistogramsOn = kFALSE;
257 for (Int_t i=0; i<2; ++i)
259 fhNClustersITS[i] = 0;
260 fhNClustersTPC[i] = 0;
262 fhChi2PerClusterITS[i] = 0;
263 fhChi2PerClusterTPC[i] = 0;
275 fhDXYNormalized[i] = 0;
276 fhDZNormalized[i] = 0;
277 fhDXYvsDZNormalized[i] = 0;
278 fhNSigmaToVertex[i] = 0;
283 fhCutCorrelation = 0;
286 //_____________________________________________________________________________
287 AliESDtrackCuts &AliESDtrackCuts::operator=(const AliESDtrackCuts &c)
290 // Assignment operator
293 if (this != &c) ((AliESDtrackCuts &) c).Copy(*this);
297 //_____________________________________________________________________________
298 void AliESDtrackCuts::Copy(TObject &c) const
304 AliESDtrackCuts& target = (AliESDtrackCuts &) c;
308 target.fCutMinNClusterTPC = fCutMinNClusterTPC;
309 target.fCutMinNClusterITS = fCutMinNClusterITS;
311 target.fCutMaxChi2PerClusterTPC = fCutMaxChi2PerClusterTPC;
312 target.fCutMaxChi2PerClusterITS = fCutMaxChi2PerClusterITS;
314 target.fCutMaxC11 = fCutMaxC11;
315 target.fCutMaxC22 = fCutMaxC22;
316 target.fCutMaxC33 = fCutMaxC33;
317 target.fCutMaxC44 = fCutMaxC44;
318 target.fCutMaxC55 = fCutMaxC55;
320 target.fCutAcceptKinkDaughters = fCutAcceptKinkDaughters;
321 target.fCutRequireTPCRefit = fCutRequireTPCRefit;
322 target.fCutRequireITSRefit = fCutRequireITSRefit;
324 target.fCutNsigmaToVertex = fCutNsigmaToVertex;
325 target.fCutSigmaToVertexRequired = fCutSigmaToVertexRequired;
327 target.fPMin = fPMin;
328 target.fPMax = fPMax;
329 target.fPtMin = fPtMin;
330 target.fPtMax = fPtMax;
331 target.fPxMin = fPxMin;
332 target.fPxMax = fPxMax;
333 target.fPyMin = fPyMin;
334 target.fPyMax = fPyMax;
335 target.fPzMin = fPzMin;
336 target.fPzMax = fPzMax;
337 target.fEtaMin = fEtaMin;
338 target.fEtaMax = fEtaMax;
339 target.fRapMin = fRapMin;
340 target.fRapMax = fRapMax;
342 target.fHistogramsOn = fHistogramsOn;
344 for (Int_t i=0; i<2; ++i)
346 if (fhNClustersITS[i]) target.fhNClustersITS[i] = (TH1F*) fhNClustersITS[i]->Clone();
347 if (fhNClustersTPC[i]) target.fhNClustersTPC[i] = (TH1F*) fhNClustersTPC[i]->Clone();
349 if (fhChi2PerClusterITS[i]) target.fhChi2PerClusterITS[i] = (TH1F*) fhChi2PerClusterITS[i]->Clone();
350 if (fhChi2PerClusterTPC[i]) target.fhChi2PerClusterTPC[i] = (TH1F*) fhChi2PerClusterTPC[i]->Clone();
352 if (fhC11[i]) target.fhC11[i] = (TH1F*) fhC11[i]->Clone();
353 if (fhC22[i]) target.fhC22[i] = (TH1F*) fhC22[i]->Clone();
354 if (fhC33[i]) target.fhC33[i] = (TH1F*) fhC33[i]->Clone();
355 if (fhC44[i]) target.fhC44[i] = (TH1F*) fhC44[i]->Clone();
356 if (fhC55[i]) target.fhC55[i] = (TH1F*) fhC55[i]->Clone();
358 if (fhDXY[i]) target.fhDXY[i] = (TH1F*) fhDXY[i]->Clone();
359 if (fhDZ[i]) target.fhDZ[i] = (TH1F*) fhDZ[i]->Clone();
360 if (fhDXYvsDZ[i]) target.fhDXYvsDZ[i] = (TH2F*) fhDXYvsDZ[i]->Clone();
362 if (fhDXYNormalized[i]) target.fhDXYNormalized[i] = (TH1F*) fhDXYNormalized[i]->Clone();
363 if (fhDZNormalized[i]) target.fhDZNormalized[i] = (TH1F*) fhDZNormalized[i]->Clone();
364 if (fhDXYvsDZNormalized[i]) target.fhDXYvsDZNormalized[i] = (TH2F*) fhDXYvsDZNormalized[i]->Clone();
365 if (fhNSigmaToVertex[i]) target.fhNSigmaToVertex[i] = (TH1F*) fhNSigmaToVertex[i]->Clone();
367 if (ffDTheoretical) target.ffDTheoretical = (TF1*) ffDTheoretical->Clone();
369 if (fhCutStatistics) target.fhCutStatistics = (TH1F*) fhCutStatistics->Clone();
370 if (fhCutCorrelation) target.fhCutCorrelation = (TH2F*) fhCutCorrelation->Clone();
375 //_____________________________________________________________________________
376 Long64_t AliESDtrackCuts::Merge(TCollection* list) {
377 // Merge a list of AliESDtrackCuts objects with this (needed for PROOF)
378 // Returns the number of merged objects (including this)
389 TIterator* iter = list->MakeIterator();
393 // collection of measured and generated histograms
395 while ((obj = iter->Next())) {
397 AliESDtrackCuts* entry = dynamic_cast<AliESDtrackCuts*>(obj);
401 if (!entry->fHistogramsOn)
404 for (Int_t i=0; i<2; i++) {
406 fhNClustersITS[i] ->Add(entry->fhNClustersITS[i] );
407 fhNClustersTPC[i] ->Add(entry->fhNClustersTPC[i] );
409 fhChi2PerClusterITS[i] ->Add(entry->fhChi2PerClusterITS[i]);
410 fhChi2PerClusterTPC[i] ->Add(entry->fhChi2PerClusterTPC[i]);
412 fhC11[i] ->Add(entry->fhC11[i] );
413 fhC22[i] ->Add(entry->fhC22[i] );
414 fhC33[i] ->Add(entry->fhC33[i] );
415 fhC44[i] ->Add(entry->fhC44[i] );
416 fhC55[i] ->Add(entry->fhC55[i] );
418 fhDXY[i] ->Add(entry->fhDXY[i] );
419 fhDZ[i] ->Add(entry->fhDZ[i] );
420 fhDXYvsDZ[i] ->Add(entry->fhDXYvsDZ[i] );
422 fhDXYNormalized[i] ->Add(entry->fhDXYNormalized[i] );
423 fhDZNormalized[i] ->Add(entry->fhDZNormalized[i] );
424 fhDXYvsDZNormalized[i] ->Add(entry->fhDXYvsDZNormalized[i]);
425 fhNSigmaToVertex[i] ->Add(entry->fhNSigmaToVertex[i]);
429 fhCutStatistics ->Add(entry->fhCutStatistics);
430 fhCutCorrelation ->Add(entry->fhCutCorrelation);
439 //____________________________________________________________________
440 Float_t AliESDtrackCuts::GetSigmaToVertex(AliESDtrack* esdTrack)
442 // Calculates the number of sigma to the vertex.
447 esdTrack->GetImpactParameters(b,bCov);
448 if (bCov[0]<=0 || bCov[2]<=0) {
449 AliDebug(1, "Estimated b resolution lower or equal zero!");
450 bCov[0]=0; bCov[2]=0;
452 bRes[0] = TMath::Sqrt(bCov[0]);
453 bRes[1] = TMath::Sqrt(bCov[2]);
455 // -----------------------------------
456 // How to get to a n-sigma cut?
458 // The accumulated statistics from 0 to d is
460 // -> Erf(d/Sqrt(2)) for a 1-dim gauss (d = n_sigma)
461 // -> 1 - Exp(-d**2) for a 2-dim gauss (d*d = dx*dx + dy*dy != n_sigma)
463 // It means that for a 2-dim gauss: n_sigma(d) = Sqrt(2)*ErfInv(1 - Exp((-x**2)/2)
464 // Can this be expressed in a different way?
466 if (bRes[0] == 0 || bRes[1] ==0)
469 Float_t d = TMath::Sqrt(TMath::Power(b[0]/bRes[0],2) + TMath::Power(b[1]/bRes[1],2));
471 // stupid rounding problem screws up everything:
472 // if d is too big, TMath::Exp(...) gets 0, and TMath::ErfInverse(1) that should be infinite, gets 0 :(
473 if (TMath::Exp(-d * d / 2) < 1e-10)
476 d = TMath::ErfInverse(1 - TMath::Exp(-d * d / 2)) * TMath::Sqrt(2);
480 void AliESDtrackCuts::EnableNeededBranches(TTree* tree)
482 // enables the branches needed by AcceptTrack, for a list see comment of AcceptTrack
484 tree->SetBranchStatus("fTracks.fFlags", 1);
485 tree->SetBranchStatus("fTracks.fITSncls", 1);
486 tree->SetBranchStatus("fTracks.fTPCncls", 1);
487 tree->SetBranchStatus("fTracks.fITSchi2", 1);
488 tree->SetBranchStatus("fTracks.fTPCchi2", 1);
489 tree->SetBranchStatus("fTracks.fC*", 1);
490 tree->SetBranchStatus("fTracks.fD", 1);
491 tree->SetBranchStatus("fTracks.fZ", 1);
492 tree->SetBranchStatus("fTracks.fCdd", 1);
493 tree->SetBranchStatus("fTracks.fCdz", 1);
494 tree->SetBranchStatus("fTracks.fCzz", 1);
495 tree->SetBranchStatus("fTracks.fP*", 1);
496 tree->SetBranchStatus("fTracks.fR*", 1);
497 tree->SetBranchStatus("fTracks.fKinkIndexes*", 1);
500 //____________________________________________________________________
502 AliESDtrackCuts::AcceptTrack(AliESDtrack* esdTrack) {
504 // figure out if the tracks survives all the track cuts defined
506 // the different quality parameter and kinematic values are first
507 // retrieved from the track. then it is found out what cuts the
508 // track did not survive and finally the cuts are imposed.
510 // this function needs the following branches:
516 // fTracks.fC //GetExternalCovariance
517 // fTracks.fD //GetImpactParameters
518 // fTracks.fZ //GetImpactParameters
519 // fTracks.fCdd //GetImpactParameters
520 // fTracks.fCdz //GetImpactParameters
521 // fTracks.fCzz //GetImpactParameters
522 // fTracks.fP //GetPxPyPz
523 // fTracks.fR //GetMass
524 // fTracks.fP //GetMass
525 // fTracks.fKinkIndexes
527 UInt_t status = esdTrack->GetStatus();
532 // getting quality parameters from the ESD track
533 Int_t nClustersITS = esdTrack->GetITSclusters(fIdxInt);
534 Int_t nClustersTPC = esdTrack->GetTPCclusters(fIdxInt);
538 Float_t chi2PerClusterITS = -1;
539 Float_t chi2PerClusterTPC = -1;
541 chi2PerClusterITS = esdTrack->GetITSchi2()/Float_t(nClustersITS);
543 chi2PerClusterTPC = esdTrack->GetTPCchi2()/Float_t(nClustersTPC);
546 esdTrack->GetExternalCovariance(extCov);
548 // getting the track to vertex parameters
549 Float_t nSigmaToVertex = GetSigmaToVertex(esdTrack);
551 // getting the kinematic variables of the track
552 // (assuming the mass is known)
554 esdTrack->GetPxPyPz(p);
555 Float_t momentum = TMath::Sqrt(TMath::Power(p[0],2) + TMath::Power(p[1],2) + TMath::Power(p[2],2));
556 Float_t pt = TMath::Sqrt(TMath::Power(p[0],2) + TMath::Power(p[1],2));
557 Float_t energy = TMath::Sqrt(TMath::Power(esdTrack->GetMass(),2) + TMath::Power(momentum,2));
560 //y-eta related calculations
563 if((momentum != TMath::Abs(p[2]))&&(momentum != 0))
564 eta = 0.5*TMath::Log((momentum + p[2])/(momentum - p[2]));
565 if((energy != TMath::Abs(p[2]))&&(momentum != 0))
566 y = 0.5*TMath::Log((energy + p[2])/(energy - p[2]));
569 //########################################################################
573 for (Int_t i=0; i<kNCuts; i++) cuts[i]=kFALSE;
575 // track quality cuts
576 if (fCutRequireTPCRefit && (status&AliESDtrack::kTPCrefit)==0)
578 if (fCutRequireITSRefit && (status&AliESDtrack::kITSrefit)==0)
580 if (nClustersTPC<fCutMinNClusterTPC)
582 if (nClustersITS<fCutMinNClusterITS)
584 if (chi2PerClusterTPC>fCutMaxChi2PerClusterTPC)
586 if (chi2PerClusterITS>fCutMaxChi2PerClusterITS)
588 if (extCov[0] > fCutMaxC11)
590 if (extCov[2] > fCutMaxC22)
592 if (extCov[5] > fCutMaxC33)
594 if (extCov[9] > fCutMaxC44)
596 if (extCov[14] > fCutMaxC55)
598 if (nSigmaToVertex > fCutNsigmaToVertex && fCutSigmaToVertexRequired)
600 // if n sigma could not be calculated
601 if (nSigmaToVertex<0 && fCutSigmaToVertexRequired)
603 if (!fCutAcceptKinkDaughters && esdTrack->GetKinkIndex(0)>0)
605 // track kinematics cut
606 if((momentum < fPMin) || (momentum > fPMax))
608 if((pt < fPtMin) || (pt > fPtMax))
610 if((p[0] < fPxMin) || (p[0] > fPxMax))
612 if((p[1] < fPyMin) || (p[1] > fPyMax))
614 if((p[2] < fPzMin) || (p[2] > fPzMax))
616 if((eta < fEtaMin) || (eta > fEtaMax))
618 if((y < fRapMin) || (y > fRapMax))
622 for (Int_t i=0; i<kNCuts; i++)
623 if (cuts[i]) cut = kTRUE;
625 //########################################################################
626 // filling histograms
628 fhCutStatistics->Fill(fhCutStatistics->GetBinCenter(fhCutStatistics->GetXaxis()->FindBin("n tracks")));
631 fhCutStatistics->Fill(fhCutStatistics->GetBinCenter(fhCutStatistics->GetXaxis()->FindBin("n cut tracks")));
633 for (Int_t i=0; i<kNCuts; i++) {
635 fhCutStatistics->Fill(fhCutStatistics->GetBinCenter(fhCutStatistics->GetXaxis()->FindBin(fgkCutNames[i])));
637 for (Int_t j=i; j<kNCuts; j++) {
638 if (cuts[i] && cuts[j]) {
639 Float_t x = fhCutCorrelation->GetXaxis()->GetBinCenter(fhCutCorrelation->GetXaxis()->FindBin(fgkCutNames[i]));
640 Float_t y = fhCutCorrelation->GetYaxis()->GetBinCenter(fhCutCorrelation->GetYaxis()->FindBin(fgkCutNames[j]));
641 fhCutCorrelation->Fill(x,y);
647 fhNClustersITS[0]->Fill(nClustersITS);
648 fhNClustersTPC[0]->Fill(nClustersTPC);
649 fhChi2PerClusterITS[0]->Fill(chi2PerClusterITS);
650 fhChi2PerClusterTPC[0]->Fill(chi2PerClusterTPC);
652 fhC11[0]->Fill(extCov[0]);
653 fhC22[0]->Fill(extCov[2]);
654 fhC33[0]->Fill(extCov[5]);
655 fhC44[0]->Fill(extCov[9]);
656 fhC55[0]->Fill(extCov[14]);
661 esdTrack->GetImpactParameters(b,bCov);
662 if (bCov[0]<=0 || bCov[2]<=0) {
663 AliDebug(1, "Estimated b resolution lower or equal zero!");
664 bCov[0]=0; bCov[2]=0;
666 bRes[0] = TMath::Sqrt(bCov[0]);
667 bRes[1] = TMath::Sqrt(bCov[2]);
670 fhDXY[0]->Fill(b[0]);
671 fhDXYvsDZ[0]->Fill(b[1],b[0]);
673 if (bRes[0]!=0 && bRes[1]!=0) {
674 fhDZNormalized[0]->Fill(b[1]/bRes[1]);
675 fhDXYNormalized[0]->Fill(b[0]/bRes[0]);
676 fhDXYvsDZNormalized[0]->Fill(b[1]/bRes[1], b[0]/bRes[0]);
677 fhNSigmaToVertex[0]->Fill(nSigmaToVertex);
681 //########################################################################
683 if (cut) return kFALSE;
685 //########################################################################
686 // filling histograms after cut
688 fhNClustersITS[1]->Fill(nClustersITS);
689 fhNClustersTPC[1]->Fill(nClustersTPC);
690 fhChi2PerClusterITS[1]->Fill(chi2PerClusterITS);
691 fhChi2PerClusterTPC[1]->Fill(chi2PerClusterTPC);
693 fhC11[1]->Fill(extCov[0]);
694 fhC22[1]->Fill(extCov[2]);
695 fhC33[1]->Fill(extCov[5]);
696 fhC44[1]->Fill(extCov[9]);
697 fhC55[1]->Fill(extCov[14]);
702 esdTrack->GetImpactParameters(b,bCov);
703 if (bCov[0]<=0 || bCov[2]<=0) {
704 AliDebug(1, "Estimated b resolution lower or equal zero!");
705 bCov[0]=0; bCov[2]=0;
707 bRes[0] = TMath::Sqrt(bCov[0]);
708 bRes[1] = TMath::Sqrt(bCov[2]);
711 fhDXY[1]->Fill(b[0]);
712 fhDXYvsDZ[1]->Fill(b[1],b[0]);
714 if (bRes[0]!=0 && bRes[1]!=0)
716 fhDZNormalized[1]->Fill(b[1]/bRes[1]);
717 fhDXYNormalized[1]->Fill(b[0]/bRes[0]);
718 fhDXYvsDZNormalized[1]->Fill(b[1]/bRes[1], b[0]/bRes[0]);
719 fhNSigmaToVertex[1]->Fill(nSigmaToVertex);
726 //____________________________________________________________________
727 TObjArray* AliESDtrackCuts::GetAcceptedTracks(AliESD* esd)
730 // returns an array of all tracks that pass the cuts
733 TObjArray* acceptedTracks = new TObjArray();
735 // loop over esd tracks
736 for (Int_t iTrack = 0; iTrack < esd->GetNumberOfTracks(); iTrack++) {
737 AliESDtrack* track = esd->GetTrack(iTrack);
739 if (AcceptTrack(track))
740 acceptedTracks->Add(track);
743 return acceptedTracks;
746 //____________________________________________________________________
747 Int_t AliESDtrackCuts::CountAcceptedTracks(AliESD* esd)
750 // returns an the number of tracks that pass the cuts
755 // loop over esd tracks
756 for (Int_t iTrack = 0; iTrack < esd->GetNumberOfTracks(); iTrack++) {
757 AliESDtrack* track = esd->GetTrack(iTrack);
759 if (AcceptTrack(track))
766 //____________________________________________________________________
767 TObjArray* AliESDtrackCuts::GetAcceptedTracks(AliESDEvent* esd)
770 // returns an array of all tracks that pass the cuts
773 TObjArray* acceptedTracks = new TObjArray();
775 // loop over esd tracks
776 for (Int_t iTrack = 0; iTrack < esd->GetNumberOfTracks(); iTrack++) {
777 AliESDtrack* track = esd->GetTrack(iTrack);
779 if (AcceptTrack(track))
780 acceptedTracks->Add(track);
783 return acceptedTracks;
786 //____________________________________________________________________
787 Int_t AliESDtrackCuts::CountAcceptedTracks(AliESDEvent* esd)
790 // returns an the number of tracks that pass the cuts
795 // loop over esd tracks
796 for (Int_t iTrack = 0; iTrack < esd->GetNumberOfTracks(); iTrack++) {
797 AliESDtrack* track = esd->GetTrack(iTrack);
799 if (AcceptTrack(track))
806 //____________________________________________________________________
807 void AliESDtrackCuts::DefineHistograms(Int_t color) {
809 // diagnostics histograms are defined
814 Bool_t oldStatus = TH1::AddDirectoryStatus();
815 TH1::AddDirectory(kFALSE);
817 //###################################################################################
818 // defining histograms
820 fhCutStatistics = new TH1F("cut_statistics","cut statistics",kNCuts+4,-0.5,kNCuts+3.5);
822 fhCutStatistics->GetXaxis()->SetBinLabel(1,"n tracks");
823 fhCutStatistics->GetXaxis()->SetBinLabel(2,"n cut tracks");
825 fhCutCorrelation = new TH2F("cut_correlation","cut correlation",kNCuts,-0.5,kNCuts-0.5,kNCuts,-0.5,kNCuts-0.5);;
827 for (Int_t i=0; i<kNCuts; i++) {
828 fhCutStatistics->GetXaxis()->SetBinLabel(i+4,fgkCutNames[i]);
829 fhCutCorrelation->GetXaxis()->SetBinLabel(i+1,fgkCutNames[i]);
830 fhCutCorrelation->GetYaxis()->SetBinLabel(i+1,fgkCutNames[i]);
833 fhCutStatistics ->SetLineColor(color);
834 fhCutCorrelation ->SetLineColor(color);
835 fhCutStatistics ->SetLineWidth(2);
836 fhCutCorrelation ->SetLineWidth(2);
839 for (Int_t i=0; i<2; i++) {
840 if (i==0) sprintf(str," ");
841 else sprintf(str,"_cut");
843 fhNClustersITS[i] = new TH1F(Form("nClustersITS%s",str) ,"",8,-0.5,7.5);
844 fhNClustersTPC[i] = new TH1F(Form("nClustersTPC%s",str) ,"",165,-0.5,164.5);
845 fhChi2PerClusterITS[i] = new TH1F(Form("chi2PerClusterITS%s",str),"",500,0,10);
846 fhChi2PerClusterTPC[i] = new TH1F(Form("chi2PerClusterTPC%s",str),"",500,0,10);
848 fhC11[i] = new TH1F(Form("covMatrixDiagonal11%s",str),"",2000,0,20);
849 fhC22[i] = new TH1F(Form("covMatrixDiagonal22%s",str),"",2000,0,20);
850 fhC33[i] = new TH1F(Form("covMatrixDiagonal33%s",str),"",1000,0,1);
851 fhC44[i] = new TH1F(Form("covMatrixDiagonal44%s",str),"",1000,0,5);
852 fhC55[i] = new TH1F(Form("covMatrixDiagonal55%s",str),"",1000,0,5);
854 fhDXY[i] = new TH1F(Form("dXY%s",str) ,"",500,-10,10);
855 fhDZ[i] = new TH1F(Form("dZ%s",str) ,"",500,-10,10);
856 fhDXYvsDZ[i] = new TH2F(Form("dXYvsDZ%s",str),"",200,-10,10,200,-10,10);
858 fhDXYNormalized[i] = new TH1F(Form("dXYNormalized%s",str) ,"",500,-10,10);
859 fhDZNormalized[i] = new TH1F(Form("dZNormalized%s",str) ,"",500,-10,10);
860 fhDXYvsDZNormalized[i] = new TH2F(Form("dXYvsDZNormalized%s",str),"",200,-10,10,200,-10,10);
862 fhNSigmaToVertex[i] = new TH1F(Form("nSigmaToVertex%s",str),"",500,0,50);
864 fhNClustersITS[i]->SetTitle("n ITS clusters");
865 fhNClustersTPC[i]->SetTitle("n TPC clusters");
866 fhChi2PerClusterITS[i]->SetTitle("#Chi^{2} per ITS cluster");
867 fhChi2PerClusterTPC[i]->SetTitle("#Chi^{2} per TPC cluster");
869 fhC11[i]->SetTitle("cov 11 : #sigma_{y}^{2} [cm^{2}]");
870 fhC22[i]->SetTitle("cov 22 : #sigma_{z}^{2} [cm^{2}]");
871 fhC33[i]->SetTitle("cov 33 : #sigma_{sin(#phi)}^{2}");
872 fhC44[i]->SetTitle("cov 44 : #sigma_{tan(#theta_{dip})}^{2}");
873 fhC55[i]->SetTitle("cov 55 : #sigma_{1/p_{T}}^{2} [(c/GeV)^2]");
875 fhDXY[i]->SetTitle("transverse impact parameter");
876 fhDZ[i]->SetTitle("longitudinal impact parameter");
877 fhDXYvsDZ[i]->SetTitle("longitudinal impact parameter");
878 fhDXYvsDZ[i]->SetYTitle("transverse impact parameter");
880 fhDXYNormalized[i]->SetTitle("normalized trans impact par");
881 fhDZNormalized[i]->SetTitle("normalized long impact par");
882 fhDXYvsDZNormalized[i]->SetTitle("normalized long impact par");
883 fhDXYvsDZNormalized[i]->SetYTitle("normalized trans impact par");
884 fhNSigmaToVertex[i]->SetTitle("n #sigma to vertex");
886 fhNClustersITS[i]->SetLineColor(color); fhNClustersITS[i]->SetLineWidth(2);
887 fhNClustersTPC[i]->SetLineColor(color); fhNClustersTPC[i]->SetLineWidth(2);
888 fhChi2PerClusterITS[i]->SetLineColor(color); fhChi2PerClusterITS[i]->SetLineWidth(2);
889 fhChi2PerClusterTPC[i]->SetLineColor(color); fhChi2PerClusterTPC[i]->SetLineWidth(2);
891 fhC11[i]->SetLineColor(color); fhC11[i]->SetLineWidth(2);
892 fhC22[i]->SetLineColor(color); fhC22[i]->SetLineWidth(2);
893 fhC33[i]->SetLineColor(color); fhC33[i]->SetLineWidth(2);
894 fhC44[i]->SetLineColor(color); fhC44[i]->SetLineWidth(2);
895 fhC55[i]->SetLineColor(color); fhC55[i]->SetLineWidth(2);
897 fhDXY[i]->SetLineColor(color); fhDXY[i]->SetLineWidth(2);
898 fhDZ[i]->SetLineColor(color); fhDZ[i]->SetLineWidth(2);
900 fhDXYNormalized[i]->SetLineColor(color); fhDXYNormalized[i]->SetLineWidth(2);
901 fhDZNormalized[i]->SetLineColor(color); fhDZNormalized[i]->SetLineWidth(2);
902 fhNSigmaToVertex[i]->SetLineColor(color); fhNSigmaToVertex[i]->SetLineWidth(2);
905 // The number of sigmas to the vertex is per definition gaussian
906 ffDTheoretical = new TF1("nSigmaToVertexTheoretical","([0]/2.506628274)*exp(-(x**2)/2)",0,50);
907 ffDTheoretical->SetParameter(0,1);
909 TH1::AddDirectory(oldStatus);
912 //____________________________________________________________________
913 Bool_t AliESDtrackCuts::LoadHistograms(const Char_t* dir)
916 // loads the histograms from a file
917 // if dir is empty a directory with the name of this object is taken (like in SaveHistogram)
923 if (!gDirectory->cd(dir))
926 ffDTheoretical = dynamic_cast<TF1*> (gDirectory->Get("nSigmaToVertexTheory"));
928 fhCutStatistics = dynamic_cast<TH1F*> (gDirectory->Get("cut_statistics"));
929 fhCutCorrelation = dynamic_cast<TH2F*> (gDirectory->Get("cut_correlation"));
932 for (Int_t i=0; i<2; i++) {
935 gDirectory->cd("before_cuts");
940 gDirectory->cd("after_cuts");
944 fhNClustersITS[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("nClustersITS%s",str) ));
945 fhNClustersTPC[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("nClustersTPC%s",str) ));
946 fhChi2PerClusterITS[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("chi2PerClusterITS%s",str)));
947 fhChi2PerClusterTPC[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("chi2PerClusterTPC%s",str)));
949 fhC11[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("covMatrixDiagonal11%s",str)));
950 fhC22[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("covMatrixDiagonal22%s",str)));
951 fhC33[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("covMatrixDiagonal33%s",str)));
952 fhC44[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("covMatrixDiagonal44%s",str)));
953 fhC55[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("covMatrixDiagonal55%s",str)));
955 fhDXY[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("dXY%s",str) ));
956 fhDZ[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("dZ%s",str) ));
957 fhDXYvsDZ[i] = dynamic_cast<TH2F*> (gDirectory->Get(Form("dXYvsDZ%s",str)));
959 fhDXYNormalized[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("dXYNormalized%s",str) ));
960 fhDZNormalized[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("dZNormalized%s",str) ));
961 fhDXYvsDZNormalized[i] = dynamic_cast<TH2F*> (gDirectory->Get(Form("dXYvsDZNormalized%s",str)));
963 fhNSigmaToVertex[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("nSigmaToVertex%s",str)));
965 // TODO only temporary
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]->SetTitle("transverse impact parameter");
978 fhDZ[i]->SetTitle("longitudinal impact parameter");
979 fhDXYvsDZ[i]->SetTitle("longitudinal impact parameter");
980 fhDXYvsDZ[i]->SetYTitle("transverse impact parameter");
982 fhDXYNormalized[i]->SetTitle("normalized trans impact par");
983 fhDZNormalized[i]->SetTitle("normalized long impact par");
984 fhDXYvsDZNormalized[i]->SetTitle("normalized long impact par");
985 fhDXYvsDZNormalized[i]->SetYTitle("normalized trans impact par");
986 fhNSigmaToVertex[i]->SetTitle("n #sigma to vertex");*/
988 gDirectory->cd("../");
991 gDirectory->cd("..");
996 //____________________________________________________________________
997 void AliESDtrackCuts::SaveHistograms(const Char_t* dir) {
999 // saves the histograms in a directory (dir)
1002 if (!fHistogramsOn) {
1003 AliDebug(0, "Histograms not on - cannot save histograms!!!");
1010 gDirectory->mkdir(dir);
1011 gDirectory->cd(dir);
1013 gDirectory->mkdir("before_cuts");
1014 gDirectory->mkdir("after_cuts");
1016 // a factor of 2 is needed since n sigma is positive
1017 ffDTheoretical->SetParameter(0,2*fhNSigmaToVertex[0]->Integral("width"));
1018 ffDTheoretical->Write("nSigmaToVertexTheory");
1020 fhCutStatistics->Write();
1021 fhCutCorrelation->Write();
1023 for (Int_t i=0; i<2; i++) {
1025 gDirectory->cd("before_cuts");
1027 gDirectory->cd("after_cuts");
1029 fhNClustersITS[i] ->Write();
1030 fhNClustersTPC[i] ->Write();
1031 fhChi2PerClusterITS[i] ->Write();
1032 fhChi2PerClusterTPC[i] ->Write();
1042 fhDXYvsDZ[i] ->Write();
1044 fhDXYNormalized[i] ->Write();
1045 fhDZNormalized[i] ->Write();
1046 fhDXYvsDZNormalized[i] ->Write();
1047 fhNSigmaToVertex[i] ->Write();
1049 gDirectory->cd("../");
1052 gDirectory->cd("../");
1055 //____________________________________________________________________
1056 void AliESDtrackCuts::DrawHistograms()
1058 // draws some histograms
1060 TCanvas* canvas1 = new TCanvas(Form("%s_1", GetName()), "Track Quality Results1", 800, 800);
1061 canvas1->Divide(2, 2);
1064 fhNClustersTPC[0]->SetStats(kFALSE);
1065 fhNClustersTPC[0]->Draw();
1068 fhChi2PerClusterTPC[0]->SetStats(kFALSE);
1069 fhChi2PerClusterTPC[0]->Draw();
1072 fhNSigmaToVertex[0]->SetStats(kFALSE);
1073 fhNSigmaToVertex[0]->GetXaxis()->SetRangeUser(0, 10);
1074 fhNSigmaToVertex[0]->Draw();
1076 canvas1->SaveAs(Form("%s_%s.gif", GetName(), canvas1->GetName()));
1078 TCanvas* canvas2 = new TCanvas(Form("%s_2", GetName()), "Track Quality Results2", 1200, 800);
1079 canvas2->Divide(3, 2);
1082 fhC11[0]->SetStats(kFALSE);
1087 fhC22[0]->SetStats(kFALSE);
1092 fhC33[0]->SetStats(kFALSE);
1097 fhC44[0]->SetStats(kFALSE);
1102 fhC55[0]->SetStats(kFALSE);
1106 canvas2->SaveAs(Form("%s_%s.gif", GetName(), canvas2->GetName()));
1108 TCanvas* canvas3 = new TCanvas(Form("%s_3", GetName()), "Track Quality Results3", 1200, 800);
1109 canvas3->Divide(3, 2);
1112 fhDXY[0]->SetStats(kFALSE);
1117 fhDZ[0]->SetStats(kFALSE);
1122 fhDXYvsDZ[0]->SetStats(kFALSE);
1124 gPad->SetRightMargin(0.15);
1125 fhDXYvsDZ[0]->Draw("COLZ");
1128 fhDXYNormalized[0]->SetStats(kFALSE);
1130 fhDXYNormalized[0]->Draw();
1133 fhDZNormalized[0]->SetStats(kFALSE);
1135 fhDZNormalized[0]->Draw();
1138 fhDXYvsDZNormalized[0]->SetStats(kFALSE);
1140 gPad->SetRightMargin(0.15);
1141 fhDXYvsDZNormalized[0]->Draw("COLZ");
1143 canvas3->SaveAs(Form("%s_%s.gif", GetName(), canvas3->GetName()));
1145 TCanvas* canvas4 = new TCanvas(Form("%s_4", GetName()), "Track Quality Results4", 800, 500);
1146 canvas4->Divide(2, 1);
1149 fhCutStatistics->SetStats(kFALSE);
1150 fhCutStatistics->LabelsOption("v");
1151 gPad->SetBottomMargin(0.3);
1152 fhCutStatistics->Draw();
1155 fhCutCorrelation->SetStats(kFALSE);
1156 fhCutCorrelation->LabelsOption("v");
1157 gPad->SetBottomMargin(0.3);
1158 gPad->SetLeftMargin(0.3);
1159 fhCutCorrelation->Draw("COLZ");
1161 canvas4->SaveAs(Form("%s_%s.gif", GetName(), canvas4->GetName()));
1164 fhDXYvsDZNormalized[0]->SetStats(kFALSE);
1165 fhDXYvsDZNormalized[0]->DrawCopy("COLZ");
1168 fhNClustersTPC[0]->SetStats(kFALSE);
1169 fhNClustersTPC[0]->DrawCopy();
1172 fhChi2PerClusterITS[0]->SetStats(kFALSE);
1173 fhChi2PerClusterITS[0]->DrawCopy();
1174 fhChi2PerClusterITS[1]->SetLineColor(2);
1175 fhChi2PerClusterITS[1]->DrawCopy("SAME");
1178 fhChi2PerClusterTPC[0]->SetStats(kFALSE);
1179 fhChi2PerClusterTPC[0]->DrawCopy();
1180 fhChi2PerClusterTPC[1]->SetLineColor(2);
1181 fhChi2PerClusterTPC[1]->DrawCopy("SAME");*/