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"
21 #include <AliESDtrack.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() : TNamed(),
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),
93 // default constructor
99 //____________________________________________________________________
100 AliESDtrackCuts::AliESDtrackCuts(Char_t* name, Char_t* title) : TNamed(name,title),
101 fCutMinNClusterTPC(0),
102 fCutMinNClusterITS(0),
103 fCutMaxChi2PerClusterTPC(0),
104 fCutMaxChi2PerClusterITS(0),
110 fCutAcceptKinkDaughters(0),
111 fCutRequireTPCRefit(0),
112 fCutRequireITSRefit(0),
113 fCutNsigmaToVertex(0),
114 fCutSigmaToVertexRequired(0),
139 //##############################################################################
140 // setting default cuts
141 SetMinNClustersTPC();
142 SetMinNClustersITS();
143 SetMaxChi2PerClusterTPC();
144 SetMaxChi2PerClusterITS();
145 SetMaxCovDiagonalElements();
146 SetRequireTPCRefit();
147 SetRequireITSRefit();
148 SetAcceptKingDaughters();
149 SetMinNsigmaToVertex();
150 SetRequireSigmaToVertex();
162 //_____________________________________________________________________________
163 AliESDtrackCuts::AliESDtrackCuts(const AliESDtrackCuts &c) : TNamed(c),
164 fCutMinNClusterTPC(0),
165 fCutMinNClusterITS(0),
166 fCutMaxChi2PerClusterTPC(0),
167 fCutMaxChi2PerClusterITS(0),
173 fCutAcceptKinkDaughters(0),
174 fCutRequireTPCRefit(0),
175 fCutRequireITSRefit(0),
176 fCutNsigmaToVertex(0),
177 fCutSigmaToVertexRequired(0),
201 ((AliESDtrackCuts &) c).Copy(*this);
204 AliESDtrackCuts::~AliESDtrackCuts()
210 for (Int_t i=0; i<2; i++) {
212 if (fhNClustersITS[i])
213 delete fhNClustersITS[i];
214 if (fhNClustersTPC[i])
215 delete fhNClustersTPC[i];
216 if (fhChi2PerClusterITS[i])
217 delete fhChi2PerClusterITS[i];
218 if (fhChi2PerClusterTPC[i])
219 delete fhChi2PerClusterTPC[i];
238 if (fhDXYNormalized[i])
239 delete fhDXYNormalized[i];
240 if (fhDZNormalized[i])
241 delete fhDZNormalized[i];
242 if (fhDXYvsDZNormalized[i])
243 delete fhDXYvsDZNormalized[i];
244 if (fhNSigmaToVertex[i])
245 delete fhNSigmaToVertex[i];
249 delete ffDTheoretical;
252 delete fhCutStatistics;
253 if (fhCutCorrelation)
254 delete fhCutCorrelation;
257 void AliESDtrackCuts::Init()
260 // sets everything to zero
263 fCutMinNClusterTPC = 0;
264 fCutMinNClusterITS = 0;
266 fCutMaxChi2PerClusterTPC = 0;
267 fCutMaxChi2PerClusterITS = 0;
275 fCutAcceptKinkDaughters = 0;
276 fCutRequireTPCRefit = 0;
277 fCutRequireITSRefit = 0;
279 fCutNsigmaToVertex = 0;
280 fCutSigmaToVertexRequired = 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;
317 fhDXYNormalized[i] = 0;
318 fhDZNormalized[i] = 0;
319 fhDXYvsDZNormalized[i] = 0;
320 fhNSigmaToVertex[i] = 0;
325 fhCutCorrelation = 0;
328 //_____________________________________________________________________________
329 AliESDtrackCuts &AliESDtrackCuts::operator=(const AliESDtrackCuts &c)
332 // Assignment operator
335 if (this != &c) ((AliESDtrackCuts &) c).Copy(*this);
339 //_____________________________________________________________________________
340 void AliESDtrackCuts::Copy(TObject &c) const
346 AliESDtrackCuts& target = (AliESDtrackCuts &) c;
350 target.fCutMinNClusterTPC = fCutMinNClusterTPC;
351 target.fCutMinNClusterITS = fCutMinNClusterITS;
353 target.fCutMaxChi2PerClusterTPC = fCutMaxChi2PerClusterTPC;
354 target.fCutMaxChi2PerClusterITS = fCutMaxChi2PerClusterITS;
356 target.fCutMaxC11 = fCutMaxC11;
357 target.fCutMaxC22 = fCutMaxC22;
358 target.fCutMaxC33 = fCutMaxC33;
359 target.fCutMaxC44 = fCutMaxC44;
360 target.fCutMaxC55 = fCutMaxC55;
362 target.fCutAcceptKinkDaughters = fCutAcceptKinkDaughters;
363 target.fCutRequireTPCRefit = fCutRequireTPCRefit;
364 target.fCutRequireITSRefit = fCutRequireITSRefit;
366 target.fCutNsigmaToVertex = fCutNsigmaToVertex;
367 target.fCutSigmaToVertexRequired = fCutSigmaToVertexRequired;
369 target.fPMin = fPMin;
370 target.fPMax = fPMax;
371 target.fPtMin = fPtMin;
372 target.fPtMax = fPtMax;
373 target.fPxMin = fPxMin;
374 target.fPxMax = fPxMax;
375 target.fPyMin = fPyMin;
376 target.fPyMax = fPyMax;
377 target.fPzMin = fPzMin;
378 target.fPzMax = fPzMax;
379 target.fEtaMin = fEtaMin;
380 target.fEtaMax = fEtaMax;
381 target.fRapMin = fRapMin;
382 target.fRapMax = fRapMax;
384 target.fHistogramsOn = fHistogramsOn;
386 for (Int_t i=0; i<2; ++i)
388 if (fhNClustersITS[i]) target.fhNClustersITS[i] = (TH1F*) fhNClustersITS[i]->Clone();
389 if (fhNClustersTPC[i]) target.fhNClustersTPC[i] = (TH1F*) fhNClustersTPC[i]->Clone();
391 if (fhChi2PerClusterITS[i]) target.fhChi2PerClusterITS[i] = (TH1F*) fhChi2PerClusterITS[i]->Clone();
392 if (fhChi2PerClusterTPC[i]) target.fhChi2PerClusterTPC[i] = (TH1F*) fhChi2PerClusterTPC[i]->Clone();
394 if (fhC11[i]) target.fhC11[i] = (TH1F*) fhC11[i]->Clone();
395 if (fhC22[i]) target.fhC22[i] = (TH1F*) fhC22[i]->Clone();
396 if (fhC33[i]) target.fhC33[i] = (TH1F*) fhC33[i]->Clone();
397 if (fhC44[i]) target.fhC44[i] = (TH1F*) fhC44[i]->Clone();
398 if (fhC55[i]) target.fhC55[i] = (TH1F*) fhC55[i]->Clone();
400 if (fhDXY[i]) target.fhDXY[i] = (TH1F*) fhDXY[i]->Clone();
401 if (fhDZ[i]) target.fhDZ[i] = (TH1F*) fhDZ[i]->Clone();
402 if (fhDXYvsDZ[i]) target.fhDXYvsDZ[i] = (TH2F*) fhDXYvsDZ[i]->Clone();
404 if (fhDXYNormalized[i]) target.fhDXYNormalized[i] = (TH1F*) fhDXYNormalized[i]->Clone();
405 if (fhDZNormalized[i]) target.fhDZNormalized[i] = (TH1F*) fhDZNormalized[i]->Clone();
406 if (fhDXYvsDZNormalized[i]) target.fhDXYvsDZNormalized[i] = (TH2F*) fhDXYvsDZNormalized[i]->Clone();
407 if (fhNSigmaToVertex[i]) target.fhNSigmaToVertex[i] = (TH1F*) fhNSigmaToVertex[i]->Clone();
409 if (ffDTheoretical) target.ffDTheoretical = (TF1*) ffDTheoretical->Clone();
411 if (fhCutStatistics) target.fhCutStatistics = (TH1F*) fhCutStatistics->Clone();
412 if (fhCutCorrelation) target.fhCutCorrelation = (TH2F*) fhCutCorrelation->Clone();
417 //_____________________________________________________________________________
418 Long64_t AliESDtrackCuts::Merge(TCollection* list) {
419 // Merge a list of AliESDtrackCuts objects with this (needed for PROOF)
420 // Returns the number of merged objects (including this)
431 TIterator* iter = list->MakeIterator();
435 // collection of measured and generated histograms
437 while ((obj = iter->Next())) {
439 AliESDtrackCuts* entry = dynamic_cast<AliESDtrackCuts*>(obj);
443 if (!entry->fHistogramsOn)
446 for (Int_t i=0; i<2; i++) {
448 fhNClustersITS[i] ->Add(entry->fhNClustersITS[i] );
449 fhNClustersTPC[i] ->Add(entry->fhNClustersTPC[i] );
451 fhChi2PerClusterITS[i] ->Add(entry->fhChi2PerClusterITS[i]);
452 fhChi2PerClusterTPC[i] ->Add(entry->fhChi2PerClusterTPC[i]);
454 fhC11[i] ->Add(entry->fhC11[i] );
455 fhC22[i] ->Add(entry->fhC22[i] );
456 fhC33[i] ->Add(entry->fhC33[i] );
457 fhC44[i] ->Add(entry->fhC44[i] );
458 fhC55[i] ->Add(entry->fhC55[i] );
460 fhDXY[i] ->Add(entry->fhDXY[i] );
461 fhDZ[i] ->Add(entry->fhDZ[i] );
462 fhDXYvsDZ[i] ->Add(entry->fhDXYvsDZ[i] );
464 fhDXYNormalized[i] ->Add(entry->fhDXYNormalized[i] );
465 fhDZNormalized[i] ->Add(entry->fhDZNormalized[i] );
466 fhDXYvsDZNormalized[i] ->Add(entry->fhDXYvsDZNormalized[i]);
467 fhNSigmaToVertex[i] ->Add(entry->fhNSigmaToVertex[i]);
471 fhCutStatistics ->Add(entry->fhCutStatistics);
472 fhCutCorrelation ->Add(entry->fhCutCorrelation);
481 //____________________________________________________________________
482 Float_t AliESDtrackCuts::GetSigmaToVertex(AliESDtrack* esdTrack)
484 // Calculates the number of sigma to the vertex.
489 esdTrack->GetImpactParameters(b,bCov);
490 if (bCov[0]<=0 || bCov[2]<=0) {
491 AliDebug(1, "Estimated b resolution lower or equal zero!");
492 bCov[0]=0; bCov[2]=0;
494 bRes[0] = TMath::Sqrt(bCov[0]);
495 bRes[1] = TMath::Sqrt(bCov[2]);
497 // -----------------------------------
498 // How to get to a n-sigma cut?
500 // The accumulated statistics from 0 to d is
502 // -> Erf(d/Sqrt(2)) for a 1-dim gauss (d = n_sigma)
503 // -> 1 - Exp(-d**2) for a 2-dim gauss (d*d = dx*dx + dy*dy != n_sigma)
505 // It means that for a 2-dim gauss: n_sigma(d) = Sqrt(2)*ErfInv(1 - Exp((-x**2)/2)
506 // Can this be expressed in a different way?
508 if (bRes[0] == 0 || bRes[1] ==0)
511 Float_t d = TMath::Sqrt(TMath::Power(b[0]/bRes[0],2) + TMath::Power(b[1]/bRes[1],2));
513 // stupid rounding problem screws up everything:
514 // if d is too big, TMath::Exp(...) gets 0, and TMath::ErfInverse(1) that should be infinite, gets 0 :(
515 if (TMath::Exp(-d * d / 2) < 1e-10)
518 d = TMath::ErfInverse(1 - TMath::Exp(-d * d / 2)) * TMath::Sqrt(2);
522 void AliESDtrackCuts::EnableNeededBranches(TTree* tree)
524 // enables the branches needed by AcceptTrack, for a list see comment of AcceptTrack
526 tree->SetBranchStatus("fTracks.fFlags", 1);
527 tree->SetBranchStatus("fTracks.fITSncls", 1);
528 tree->SetBranchStatus("fTracks.fTPCncls", 1);
529 tree->SetBranchStatus("fTracks.fITSchi2", 1);
530 tree->SetBranchStatus("fTracks.fTPCchi2", 1);
531 tree->SetBranchStatus("fTracks.fC*", 1);
532 tree->SetBranchStatus("fTracks.fD", 1);
533 tree->SetBranchStatus("fTracks.fZ", 1);
534 tree->SetBranchStatus("fTracks.fCdd", 1);
535 tree->SetBranchStatus("fTracks.fCdz", 1);
536 tree->SetBranchStatus("fTracks.fCzz", 1);
537 tree->SetBranchStatus("fTracks.fP*", 1);
538 tree->SetBranchStatus("fTracks.fR*", 1);
539 tree->SetBranchStatus("fTracks.fKinkIndexes*", 1);
542 //____________________________________________________________________
544 AliESDtrackCuts::AcceptTrack(AliESDtrack* esdTrack) {
546 // figure out if the tracks survives all the track cuts defined
548 // the different quality parameter and kinematic values are first
549 // retrieved from the track. then it is found out what cuts the
550 // track did not survive and finally the cuts are imposed.
552 // this function needs the following branches:
558 // fTracks.fC //GetExternalCovariance
559 // fTracks.fD //GetImpactParameters
560 // fTracks.fZ //GetImpactParameters
561 // fTracks.fCdd //GetImpactParameters
562 // fTracks.fCdz //GetImpactParameters
563 // fTracks.fCzz //GetImpactParameters
564 // fTracks.fP //GetPxPyPz
565 // fTracks.fR //GetMass
566 // fTracks.fP //GetMass
567 // fTracks.fKinkIndexes
569 UInt_t status = esdTrack->GetStatus();
574 // getting quality parameters from the ESD track
575 Int_t nClustersITS = esdTrack->GetITSclusters(fIdxInt);
576 Int_t nClustersTPC = esdTrack->GetTPCclusters(fIdxInt);
580 Float_t chi2PerClusterITS = -1;
581 Float_t chi2PerClusterTPC = -1;
583 chi2PerClusterITS = esdTrack->GetITSchi2()/Float_t(nClustersITS);
585 chi2PerClusterTPC = esdTrack->GetTPCchi2()/Float_t(nClustersTPC);
588 esdTrack->GetExternalCovariance(extCov);
590 // getting the track to vertex parameters
591 Float_t nSigmaToVertex = GetSigmaToVertex(esdTrack);
593 // getting the kinematic variables of the track
594 // (assuming the mass is known)
596 esdTrack->GetPxPyPz(p);
597 Float_t momentum = TMath::Sqrt(TMath::Power(p[0],2) + TMath::Power(p[1],2) + TMath::Power(p[2],2));
598 Float_t pt = TMath::Sqrt(TMath::Power(p[0],2) + TMath::Power(p[1],2));
599 Float_t energy = TMath::Sqrt(TMath::Power(esdTrack->GetMass(),2) + TMath::Power(momentum,2));
602 //y-eta related calculations
605 if((momentum != TMath::Abs(p[2]))&&(momentum != 0))
606 eta = 0.5*TMath::Log((momentum + p[2])/(momentum - p[2]));
607 if((energy != TMath::Abs(p[2]))&&(momentum != 0))
608 y = 0.5*TMath::Log((energy + p[2])/(energy - p[2]));
611 //########################################################################
615 for (Int_t i=0; i<kNCuts; i++) cuts[i]=kFALSE;
617 // track quality cuts
618 if (fCutRequireTPCRefit && (status&AliESDtrack::kTPCrefit)==0)
620 if (fCutRequireITSRefit && (status&AliESDtrack::kITSrefit)==0)
622 if (nClustersTPC<fCutMinNClusterTPC)
624 if (nClustersITS<fCutMinNClusterITS)
626 if (chi2PerClusterTPC>fCutMaxChi2PerClusterTPC)
628 if (chi2PerClusterITS>fCutMaxChi2PerClusterITS)
630 if (extCov[0] > fCutMaxC11)
632 if (extCov[2] > fCutMaxC22)
634 if (extCov[5] > fCutMaxC33)
636 if (extCov[9] > fCutMaxC44)
638 if (extCov[14] > fCutMaxC55)
640 if (nSigmaToVertex > fCutNsigmaToVertex)
642 // if n sigma could not be calculated
643 if (nSigmaToVertex<0 && fCutSigmaToVertexRequired)
645 if (!fCutAcceptKinkDaughters && esdTrack->GetKinkIndex(0)>0)
647 // track kinematics cut
648 if((momentum < fPMin) || (momentum > fPMax))
650 if((pt < fPtMin) || (pt > fPtMax))
652 if((p[0] < fPxMin) || (p[0] > fPxMax))
654 if((p[1] < fPyMin) || (p[1] > fPyMax))
656 if((p[2] < fPzMin) || (p[2] > fPzMax))
658 if((eta < fEtaMin) || (eta > fEtaMax))
660 if((y < fRapMin) || (y > fRapMax))
664 for (Int_t i=0; i<kNCuts; i++)
665 if (cuts[i]) cut = kTRUE;
667 //########################################################################
668 // filling histograms
670 fhCutStatistics->Fill(fhCutStatistics->GetBinCenter(fhCutStatistics->GetXaxis()->FindBin("n tracks")));
673 fhCutStatistics->Fill(fhCutStatistics->GetBinCenter(fhCutStatistics->GetXaxis()->FindBin("n cut tracks")));
675 for (Int_t i=0; i<kNCuts; i++) {
677 fhCutStatistics->Fill(fhCutStatistics->GetBinCenter(fhCutStatistics->GetXaxis()->FindBin(fgkCutNames[i])));
679 for (Int_t j=i; j<kNCuts; j++) {
680 if (cuts[i] && cuts[j]) {
681 Float_t x = fhCutCorrelation->GetXaxis()->GetBinCenter(fhCutCorrelation->GetXaxis()->FindBin(fgkCutNames[i]));
682 Float_t y = fhCutCorrelation->GetYaxis()->GetBinCenter(fhCutCorrelation->GetYaxis()->FindBin(fgkCutNames[j]));
683 fhCutCorrelation->Fill(x,y);
689 fhNClustersITS[0]->Fill(nClustersITS);
690 fhNClustersTPC[0]->Fill(nClustersTPC);
691 fhChi2PerClusterITS[0]->Fill(chi2PerClusterITS);
692 fhChi2PerClusterTPC[0]->Fill(chi2PerClusterTPC);
694 fhC11[0]->Fill(extCov[0]);
695 fhC22[0]->Fill(extCov[2]);
696 fhC33[0]->Fill(extCov[5]);
697 fhC44[0]->Fill(extCov[9]);
698 fhC55[0]->Fill(extCov[14]);
703 esdTrack->GetImpactParameters(b,bCov);
704 if (bCov[0]<=0 || bCov[2]<=0) {
705 AliDebug(1, "Estimated b resolution lower or equal zero!");
706 bCov[0]=0; bCov[2]=0;
708 bRes[0] = TMath::Sqrt(bCov[0]);
709 bRes[1] = TMath::Sqrt(bCov[2]);
712 fhDXY[0]->Fill(b[0]);
713 fhDXYvsDZ[0]->Fill(b[1],b[0]);
715 if (bRes[0]!=0 && bRes[1]!=0) {
716 fhDZNormalized[0]->Fill(b[1]/bRes[1]);
717 fhDXYNormalized[0]->Fill(b[0]/bRes[0]);
718 fhDXYvsDZNormalized[0]->Fill(b[1]/bRes[1], b[0]/bRes[0]);
719 fhNSigmaToVertex[0]->Fill(nSigmaToVertex);
723 //########################################################################
725 if (cut) return kFALSE;
727 //########################################################################
728 // filling histograms after cut
730 fhNClustersITS[1]->Fill(nClustersITS);
731 fhNClustersTPC[1]->Fill(nClustersTPC);
732 fhChi2PerClusterITS[1]->Fill(chi2PerClusterITS);
733 fhChi2PerClusterTPC[1]->Fill(chi2PerClusterTPC);
735 fhC11[1]->Fill(extCov[0]);
736 fhC22[1]->Fill(extCov[2]);
737 fhC33[1]->Fill(extCov[5]);
738 fhC44[1]->Fill(extCov[9]);
739 fhC55[1]->Fill(extCov[14]);
744 esdTrack->GetImpactParameters(b,bCov);
745 if (bCov[0]<=0 || bCov[2]<=0) {
746 AliDebug(1, "Estimated b resolution lower or equal zero!");
747 bCov[0]=0; bCov[2]=0;
749 bRes[0] = TMath::Sqrt(bCov[0]);
750 bRes[1] = TMath::Sqrt(bCov[2]);
753 fhDXY[1]->Fill(b[0]);
754 fhDXYvsDZ[1]->Fill(b[1],b[0]);
756 if (bRes[0]!=0 && bRes[1]!=0)
758 fhDZNormalized[1]->Fill(b[1]/bRes[1]);
759 fhDXYNormalized[1]->Fill(b[0]/bRes[0]);
760 fhDXYvsDZNormalized[1]->Fill(b[1]/bRes[1], b[0]/bRes[0]);
761 fhNSigmaToVertex[1]->Fill(nSigmaToVertex);
768 //____________________________________________________________________
770 AliESDtrackCuts::GetAcceptedTracks(AliESD* esd)
773 // returns an array of all tracks that pass the cuts
776 TObjArray* acceptedTracks = new TObjArray();
778 // loop over esd tracks
779 for (Int_t iTrack = 0; iTrack < esd->GetNumberOfTracks(); iTrack++) {
780 AliESDtrack* track = esd->GetTrack(iTrack);
782 if (AcceptTrack(track))
783 acceptedTracks->Add(track);
786 return acceptedTracks;
789 //____________________________________________________________________
791 AliESDtrackCuts::CountAcceptedTracks(AliESD* esd)
794 // returns an the number of tracks that pass the cuts
799 // loop over esd tracks
800 for (Int_t iTrack = 0; iTrack < esd->GetNumberOfTracks(); iTrack++) {
801 AliESDtrack* track = esd->GetTrack(iTrack);
803 if (AcceptTrack(track))
810 //____________________________________________________________________
811 void AliESDtrackCuts::DefineHistograms(Int_t color) {
813 // diagnostics histograms are defined
818 Bool_t oldStatus = TH1::AddDirectoryStatus();
819 TH1::AddDirectory(kFALSE);
821 //###################################################################################
822 // defining histograms
824 fhCutStatistics = new TH1F("cut_statistics","cut statistics",kNCuts+4,-0.5,kNCuts+3.5);
826 fhCutStatistics->GetXaxis()->SetBinLabel(1,"n tracks");
827 fhCutStatistics->GetXaxis()->SetBinLabel(2,"n cut tracks");
829 fhCutCorrelation = new TH2F("cut_correlation","cut correlation",kNCuts,-0.5,kNCuts-0.5,kNCuts,-0.5,kNCuts-0.5);;
831 for (Int_t i=0; i<kNCuts; i++) {
832 fhCutStatistics->GetXaxis()->SetBinLabel(i+4,fgkCutNames[i]);
833 fhCutCorrelation->GetXaxis()->SetBinLabel(i+1,fgkCutNames[i]);
834 fhCutCorrelation->GetYaxis()->SetBinLabel(i+1,fgkCutNames[i]);
837 fhCutStatistics ->SetLineColor(color);
838 fhCutCorrelation ->SetLineColor(color);
839 fhCutStatistics ->SetLineWidth(2);
840 fhCutCorrelation ->SetLineWidth(2);
843 for (Int_t i=0; i<2; i++) {
844 if (i==0) sprintf(str," ");
845 else sprintf(str,"_cut");
847 fhNClustersITS[i] = new TH1F(Form("nClustersITS%s",str) ,"",8,-0.5,7.5);
848 fhNClustersTPC[i] = new TH1F(Form("nClustersTPC%s",str) ,"",165,-0.5,164.5);
849 fhChi2PerClusterITS[i] = new TH1F(Form("chi2PerClusterITS%s",str),"",500,0,10);
850 fhChi2PerClusterTPC[i] = new TH1F(Form("chi2PerClusterTPC%s",str),"",500,0,10);
852 fhC11[i] = new TH1F(Form("covMatrixDiagonal11%s",str),"",2000,0,20);
853 fhC22[i] = new TH1F(Form("covMatrixDiagonal22%s",str),"",2000,0,20);
854 fhC33[i] = new TH1F(Form("covMatrixDiagonal33%s",str),"",1000,0,1);
855 fhC44[i] = new TH1F(Form("covMatrixDiagonal44%s",str),"",1000,0,5);
856 fhC55[i] = new TH1F(Form("covMatrixDiagonal55%s",str),"",1000,0,5);
858 fhDXY[i] = new TH1F(Form("dXY%s",str) ,"",500,-10,10);
859 fhDZ[i] = new TH1F(Form("dZ%s",str) ,"",500,-10,10);
860 fhDXYvsDZ[i] = new TH2F(Form("dXYvsDZ%s",str),"",200,-10,10,200,-10,10);
862 fhDXYNormalized[i] = new TH1F(Form("dXYNormalized%s",str) ,"",500,-10,10);
863 fhDZNormalized[i] = new TH1F(Form("dZNormalized%s",str) ,"",500,-10,10);
864 fhDXYvsDZNormalized[i] = new TH2F(Form("dXYvsDZNormalized%s",str),"",200,-10,10,200,-10,10);
866 fhNSigmaToVertex[i] = new TH1F(Form("nSigmaToVertex%s",str),"",500,0,50);
868 fhNClustersITS[i]->SetTitle("n ITS clusters");
869 fhNClustersTPC[i]->SetTitle("n TPC clusters");
870 fhChi2PerClusterITS[i]->SetTitle("#Chi^{2} per ITS cluster");
871 fhChi2PerClusterTPC[i]->SetTitle("#Chi^{2} per TPC cluster");
873 fhC11[i]->SetTitle("cov 11 : #sigma_{y}^{2} [cm^{2}]");
874 fhC22[i]->SetTitle("cov 22 : #sigma_{z}^{2} [cm^{2}]");
875 fhC33[i]->SetTitle("cov 33 : #sigma_{sin(#phi)}^{2}");
876 fhC44[i]->SetTitle("cov 44 : #sigma_{tan(#theta_{dip})}^{2}");
877 fhC55[i]->SetTitle("cov 55 : #sigma_{1/p_{T}}^{2} [(c/GeV)^2]");
879 fhDXY[i]->SetTitle("transverse impact parameter");
880 fhDZ[i]->SetTitle("longitudinal impact parameter");
881 fhDXYvsDZ[i]->SetTitle("longitudinal impact parameter");
882 fhDXYvsDZ[i]->SetYTitle("transverse impact parameter");
884 fhDXYNormalized[i]->SetTitle("normalized trans impact par");
885 fhDZNormalized[i]->SetTitle("normalized long impact par");
886 fhDXYvsDZNormalized[i]->SetTitle("normalized long impact par");
887 fhDXYvsDZNormalized[i]->SetYTitle("normalized trans impact par");
888 fhNSigmaToVertex[i]->SetTitle("n #sigma to vertex");
890 fhNClustersITS[i]->SetLineColor(color); fhNClustersITS[i]->SetLineWidth(2);
891 fhNClustersTPC[i]->SetLineColor(color); fhNClustersTPC[i]->SetLineWidth(2);
892 fhChi2PerClusterITS[i]->SetLineColor(color); fhChi2PerClusterITS[i]->SetLineWidth(2);
893 fhChi2PerClusterTPC[i]->SetLineColor(color); fhChi2PerClusterTPC[i]->SetLineWidth(2);
895 fhC11[i]->SetLineColor(color); fhC11[i]->SetLineWidth(2);
896 fhC22[i]->SetLineColor(color); fhC22[i]->SetLineWidth(2);
897 fhC33[i]->SetLineColor(color); fhC33[i]->SetLineWidth(2);
898 fhC44[i]->SetLineColor(color); fhC44[i]->SetLineWidth(2);
899 fhC55[i]->SetLineColor(color); fhC55[i]->SetLineWidth(2);
901 fhDXY[i]->SetLineColor(color); fhDXY[i]->SetLineWidth(2);
902 fhDZ[i]->SetLineColor(color); fhDZ[i]->SetLineWidth(2);
904 fhDXYNormalized[i]->SetLineColor(color); fhDXYNormalized[i]->SetLineWidth(2);
905 fhDZNormalized[i]->SetLineColor(color); fhDZNormalized[i]->SetLineWidth(2);
906 fhNSigmaToVertex[i]->SetLineColor(color); fhNSigmaToVertex[i]->SetLineWidth(2);
909 // The number of sigmas to the vertex is per definition gaussian
910 ffDTheoretical = new TF1("nSigmaToVertexTheoretical","([0]/2.506628274)*exp(-(x**2)/2)",0,50);
911 ffDTheoretical->SetParameter(0,1);
913 TH1::AddDirectory(oldStatus);
916 //____________________________________________________________________
917 Bool_t AliESDtrackCuts::LoadHistograms(const Char_t* dir)
920 // loads the histograms from a file
921 // if dir is empty a directory with the name of this object is taken (like in SaveHistogram)
927 if (!gDirectory->cd(dir))
930 ffDTheoretical = dynamic_cast<TF1*> (gDirectory->Get("nSigmaToVertexTheory"));
932 fhCutStatistics = dynamic_cast<TH1F*> (gDirectory->Get("cut_statistics"));
933 fhCutCorrelation = dynamic_cast<TH2F*> (gDirectory->Get("cut_correlation"));
936 for (Int_t i=0; i<2; i++) {
939 gDirectory->cd("before_cuts");
944 gDirectory->cd("after_cuts");
948 fhNClustersITS[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("nClustersITS%s",str) ));
949 fhNClustersTPC[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("nClustersTPC%s",str) ));
950 fhChi2PerClusterITS[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("chi2PerClusterITS%s",str)));
951 fhChi2PerClusterTPC[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("chi2PerClusterTPC%s",str)));
953 fhC11[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("covMatrixDiagonal11%s",str)));
954 fhC22[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("covMatrixDiagonal22%s",str)));
955 fhC33[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("covMatrixDiagonal33%s",str)));
956 fhC44[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("covMatrixDiagonal44%s",str)));
957 fhC55[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("covMatrixDiagonal55%s",str)));
959 fhDXY[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("dXY%s",str) ));
960 fhDZ[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("dZ%s",str) ));
961 fhDXYvsDZ[i] = dynamic_cast<TH2F*> (gDirectory->Get(Form("dXYvsDZ%s",str)));
963 fhDXYNormalized[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("dXYNormalized%s",str) ));
964 fhDZNormalized[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("dZNormalized%s",str) ));
965 fhDXYvsDZNormalized[i] = dynamic_cast<TH2F*> (gDirectory->Get(Form("dXYvsDZNormalized%s",str)));
967 fhNSigmaToVertex[i] = dynamic_cast<TH1F*> (gDirectory->Get(Form("nSigmaToVertex%s",str)));
969 // TODO only temporary
970 /*fhNClustersITS[i]->SetTitle("n ITS clusters");
971 fhNClustersTPC[i]->SetTitle("n TPC clusters");
972 fhChi2PerClusterITS[i]->SetTitle("#Chi^{2} per ITS cluster");
973 fhChi2PerClusterTPC[i]->SetTitle("#Chi^{2} per TPC cluster");
975 fhC11[i]->SetTitle("cov 11 : #sigma_{y}^{2} [cm^{2}]");
976 fhC22[i]->SetTitle("cov 22 : #sigma_{z}^{2} [cm^{2}]");
977 fhC33[i]->SetTitle("cov 33 : #sigma_{sin(#phi)}^{2}");
978 fhC44[i]->SetTitle("cov 44 : #sigma_{tan(#theta_{dip})}^{2}");
979 fhC55[i]->SetTitle("cov 55 : #sigma_{1/p_{T}}^{2} [(c/GeV)^2]");
981 fhDXY[i]->SetTitle("transverse impact parameter");
982 fhDZ[i]->SetTitle("longitudinal impact parameter");
983 fhDXYvsDZ[i]->SetTitle("longitudinal impact parameter");
984 fhDXYvsDZ[i]->SetYTitle("transverse impact parameter");
986 fhDXYNormalized[i]->SetTitle("normalized trans impact par");
987 fhDZNormalized[i]->SetTitle("normalized long impact par");
988 fhDXYvsDZNormalized[i]->SetTitle("normalized long impact par");
989 fhDXYvsDZNormalized[i]->SetYTitle("normalized trans impact par");
990 fhNSigmaToVertex[i]->SetTitle("n #sigma to vertex");*/
992 gDirectory->cd("../");
995 gDirectory->cd("..");
1000 //____________________________________________________________________
1001 void AliESDtrackCuts::SaveHistograms(const Char_t* dir) {
1003 // saves the histograms in a directory (dir)
1006 if (!fHistogramsOn) {
1007 AliDebug(0, "Histograms not on - cannot save histograms!!!");
1014 gDirectory->mkdir(dir);
1015 gDirectory->cd(dir);
1017 gDirectory->mkdir("before_cuts");
1018 gDirectory->mkdir("after_cuts");
1020 // a factor of 2 is needed since n sigma is positive
1021 ffDTheoretical->SetParameter(0,2*fhNSigmaToVertex[0]->Integral("width"));
1022 ffDTheoretical->Write("nSigmaToVertexTheory");
1024 fhCutStatistics->Write();
1025 fhCutCorrelation->Write();
1027 for (Int_t i=0; i<2; i++) {
1029 gDirectory->cd("before_cuts");
1031 gDirectory->cd("after_cuts");
1033 fhNClustersITS[i] ->Write();
1034 fhNClustersTPC[i] ->Write();
1035 fhChi2PerClusterITS[i] ->Write();
1036 fhChi2PerClusterTPC[i] ->Write();
1046 fhDXYvsDZ[i] ->Write();
1048 fhDXYNormalized[i] ->Write();
1049 fhDZNormalized[i] ->Write();
1050 fhDXYvsDZNormalized[i] ->Write();
1051 fhNSigmaToVertex[i] ->Write();
1053 gDirectory->cd("../");
1056 gDirectory->cd("../");
1059 //____________________________________________________________________
1060 void AliESDtrackCuts::DrawHistograms()
1062 // draws some histograms
1064 TCanvas* canvas1 = new TCanvas(Form("%s_1", GetName()), "Track Quality Results1", 800, 800);
1065 canvas1->Divide(2, 2);
1068 fhNClustersTPC[0]->SetStats(kFALSE);
1069 fhNClustersTPC[0]->Draw();
1072 fhChi2PerClusterTPC[0]->SetStats(kFALSE);
1073 fhChi2PerClusterTPC[0]->Draw();
1076 fhNSigmaToVertex[0]->SetStats(kFALSE);
1077 fhNSigmaToVertex[0]->GetXaxis()->SetRangeUser(0, 10);
1078 fhNSigmaToVertex[0]->Draw();
1080 canvas1->SaveAs(Form("%s_%s.gif", GetName(), canvas1->GetName()));
1082 TCanvas* canvas2 = new TCanvas(Form("%s_2", GetName()), "Track Quality Results2", 1200, 800);
1083 canvas2->Divide(3, 2);
1086 fhC11[0]->SetStats(kFALSE);
1091 fhC22[0]->SetStats(kFALSE);
1096 fhC33[0]->SetStats(kFALSE);
1101 fhC44[0]->SetStats(kFALSE);
1106 fhC55[0]->SetStats(kFALSE);
1110 canvas2->SaveAs(Form("%s_%s.gif", GetName(), canvas2->GetName()));
1112 TCanvas* canvas3 = new TCanvas(Form("%s_3", GetName()), "Track Quality Results3", 1200, 800);
1113 canvas3->Divide(3, 2);
1116 fhDXY[0]->SetStats(kFALSE);
1121 fhDZ[0]->SetStats(kFALSE);
1126 fhDXYvsDZ[0]->SetStats(kFALSE);
1128 gPad->SetRightMargin(0.15);
1129 fhDXYvsDZ[0]->Draw("COLZ");
1132 fhDXYNormalized[0]->SetStats(kFALSE);
1134 fhDXYNormalized[0]->Draw();
1137 fhDZNormalized[0]->SetStats(kFALSE);
1139 fhDZNormalized[0]->Draw();
1142 fhDXYvsDZNormalized[0]->SetStats(kFALSE);
1144 gPad->SetRightMargin(0.15);
1145 fhDXYvsDZNormalized[0]->Draw("COLZ");
1147 canvas3->SaveAs(Form("%s_%s.gif", GetName(), canvas3->GetName()));
1149 TCanvas* canvas4 = new TCanvas(Form("%s_4", GetName()), "Track Quality Results4", 800, 500);
1150 canvas4->Divide(2, 1);
1153 fhCutStatistics->SetStats(kFALSE);
1154 fhCutStatistics->LabelsOption("v");
1155 gPad->SetBottomMargin(0.3);
1156 fhCutStatistics->Draw();
1159 fhCutCorrelation->SetStats(kFALSE);
1160 fhCutCorrelation->LabelsOption("v");
1161 gPad->SetBottomMargin(0.3);
1162 gPad->SetLeftMargin(0.3);
1163 fhCutCorrelation->Draw("COLZ");
1165 canvas4->SaveAs(Form("%s_%s.gif", GetName(), canvas4->GetName()));
1168 fhDXYvsDZNormalized[0]->SetStats(kFALSE);
1169 fhDXYvsDZNormalized[0]->DrawCopy("COLZ");
1172 fhNClustersTPC[0]->SetStats(kFALSE);
1173 fhNClustersTPC[0]->DrawCopy();
1176 fhChi2PerClusterITS[0]->SetStats(kFALSE);
1177 fhChi2PerClusterITS[0]->DrawCopy();
1178 fhChi2PerClusterITS[1]->SetLineColor(2);
1179 fhChi2PerClusterITS[1]->DrawCopy("SAME");
1182 fhChi2PerClusterTPC[0]->SetStats(kFALSE);
1183 fhChi2PerClusterTPC[0]->DrawCopy();
1184 fhChi2PerClusterTPC[1]->SetLineColor(2);
1185 fhChi2PerClusterTPC[1]->DrawCopy("SAME");*/