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>
31 //____________________________________________________________________
32 ClassImp(AliESDtrackCuts)
35 const Char_t* AliESDtrackCuts::fgkCutNames[kNCuts] = {
37 "require TPC standalone",
41 "#Chi^{2}/cluster TPC",
42 "#Chi^{2}/cluster ITS",
58 "trk-to-vtx max dca 2D absolute",
59 "trk-to-vtx max dca xy absolute",
60 "trk-to-vtx max dca z absolute",
61 "trk-to-vtx min dca 2D absolute",
62 "trk-to-vtx min dca xy absolute",
63 "trk-to-vtx min dca z absolute",
64 "SPD cluster requirement",
65 "SDD cluster requirement",
66 "SSD cluster requirement",
67 "require ITS stand-alone",
68 "rel 1/pt uncertainty",
72 //____________________________________________________________________
73 AliESDtrackCuts::AliESDtrackCuts(const Char_t* name, const Char_t* title) : AliAnalysisCuts(name,title),
74 fCutMinNClusterTPC(0),
75 fCutMinNClusterITS(0),
76 fCutMaxChi2PerClusterTPC(0),
77 fCutMaxChi2PerClusterITS(0),
83 fCutMaxRel1PtUncertainty(0),
84 fCutAcceptKinkDaughters(0),
85 fCutAcceptSharedTPCClusters(0),
86 fCutMaxFractionSharedTPCClusters(0),
87 fCutRequireTPCRefit(0),
88 fCutRequireTPCStandAlone(0),
89 fCutRequireITSRefit(0),
91 fCutRequireITSStandAlone(0),
92 fCutRequireITSpureSA(0),
93 fCutNsigmaToVertex(0),
94 fCutSigmaToVertexRequired(0),
95 fCutMaxDCAToVertexXY(0),
96 fCutMaxDCAToVertexZ(0),
97 fCutMinDCAToVertexXY(0),
98 fCutMinDCAToVertexZ(0),
99 fCutMaxDCAToVertexXYPtDep(""),
100 fCutMaxDCAToVertexZPtDep(""),
101 fCutMinDCAToVertexXYPtDep(""),
102 fCutMinDCAToVertexZPtDep(""),
103 f1CutMaxDCAToVertexXYPtDep(0x0),
104 f1CutMaxDCAToVertexZPtDep(0x0),
105 f1CutMinDCAToVertexXYPtDep(0x0),
106 f1CutMinDCAToVertexZPtDep(0x0),
107 fCutDCAToVertex2D(0),
133 //##############################################################################
134 // setting default cuts
135 SetMinNClustersTPC();
136 SetMinNClustersITS();
137 SetMaxChi2PerClusterTPC();
138 SetMaxChi2PerClusterITS();
139 SetMaxCovDiagonalElements();
140 SetMaxRel1PtUncertainty();
141 SetRequireTPCRefit();
142 SetRequireTPCStandAlone();
143 SetRequireITSRefit();
144 SetRequireITSPid(kFALSE);
145 SetRequireITSStandAlone(kFALSE);
146 SetRequireITSPureStandAlone(kFALSE);
147 SetAcceptKinkDaughters();
148 SetMaxNsigmaToVertex();
149 SetMaxDCAToVertexXY();
150 SetMaxDCAToVertexZ();
152 SetMinDCAToVertexXY();
153 SetMinDCAToVertexZ();
161 SetClusterRequirementITS(kSPD);
162 SetClusterRequirementITS(kSDD);
163 SetClusterRequirementITS(kSSD);
168 //_____________________________________________________________________________
169 AliESDtrackCuts::AliESDtrackCuts(const AliESDtrackCuts &c) : AliAnalysisCuts(c),
170 fCutMinNClusterTPC(0),
171 fCutMinNClusterITS(0),
172 fCutMaxChi2PerClusterTPC(0),
173 fCutMaxChi2PerClusterITS(0),
179 fCutMaxRel1PtUncertainty(0),
180 fCutAcceptKinkDaughters(0),
181 fCutAcceptSharedTPCClusters(0),
182 fCutMaxFractionSharedTPCClusters(0),
183 fCutRequireTPCRefit(0),
184 fCutRequireTPCStandAlone(0),
185 fCutRequireITSRefit(0),
186 fCutRequireITSPid(0),
187 fCutRequireITSStandAlone(0),
188 fCutRequireITSpureSA(0),
189 fCutNsigmaToVertex(0),
190 fCutSigmaToVertexRequired(0),
191 fCutMaxDCAToVertexXY(0),
192 fCutMaxDCAToVertexZ(0),
193 fCutMinDCAToVertexXY(0),
194 fCutMinDCAToVertexZ(0),
195 fCutMaxDCAToVertexXYPtDep(""),
196 fCutMaxDCAToVertexZPtDep(""),
197 fCutMinDCAToVertexXYPtDep(""),
198 fCutMinDCAToVertexZPtDep(""),
199 f1CutMaxDCAToVertexXYPtDep(0x0),
200 f1CutMaxDCAToVertexZPtDep(0x0),
201 f1CutMinDCAToVertexXYPtDep(0x0),
202 f1CutMinDCAToVertexZPtDep(0x0),
203 fCutDCAToVertex2D(0),
227 ((AliESDtrackCuts &) c).Copy(*this);
230 AliESDtrackCuts::~AliESDtrackCuts()
236 for (Int_t i=0; i<2; i++) {
238 if (fhNClustersITS[i])
239 delete fhNClustersITS[i];
240 if (fhNClustersTPC[i])
241 delete fhNClustersTPC[i];
242 if (fhChi2PerClusterITS[i])
243 delete fhChi2PerClusterITS[i];
244 if (fhChi2PerClusterTPC[i])
245 delete fhChi2PerClusterTPC[i];
257 if (fhRel1PtUncertainty[i])
258 delete fhRel1PtUncertainty[i];
269 if (fhDXYNormalized[i])
270 delete fhDXYNormalized[i];
271 if (fhDZNormalized[i])
272 delete fhDZNormalized[i];
273 if (fhDXYvsDZNormalized[i])
274 delete fhDXYvsDZNormalized[i];
275 if (fhNSigmaToVertex[i])
276 delete fhNSigmaToVertex[i];
283 if(f1CutMaxDCAToVertexXYPtDep)delete f1CutMaxDCAToVertexXYPtDep;
284 f1CutMaxDCAToVertexXYPtDep = 0;
285 if( f1CutMaxDCAToVertexZPtDep) delete f1CutMaxDCAToVertexZPtDep;
286 f1CutMaxDCAToVertexZPtDep = 0;
287 if( f1CutMinDCAToVertexXYPtDep)delete f1CutMinDCAToVertexXYPtDep;
288 f1CutMinDCAToVertexXYPtDep = 0;
289 if(f1CutMinDCAToVertexZPtDep)delete f1CutMinDCAToVertexZPtDep;
290 f1CutMinDCAToVertexZPtDep = 0;
294 delete ffDTheoretical;
297 delete fhCutStatistics;
298 if (fhCutCorrelation)
299 delete fhCutCorrelation;
302 void AliESDtrackCuts::Init()
305 // sets everything to zero
308 fCutMinNClusterTPC = 0;
309 fCutMinNClusterITS = 0;
311 fCutMaxChi2PerClusterTPC = 0;
312 fCutMaxChi2PerClusterITS = 0;
314 for (Int_t i = 0; i < 3; i++)
315 fCutClusterRequirementITS[i] = kOff;
323 fCutMaxRel1PtUncertainty = 0;
325 fCutAcceptKinkDaughters = 0;
326 fCutAcceptSharedTPCClusters = 0;
327 fCutMaxFractionSharedTPCClusters = 0;
328 fCutRequireTPCRefit = 0;
329 fCutRequireTPCStandAlone = 0;
330 fCutRequireITSRefit = 0;
331 fCutRequireITSPid = 0;
332 fCutRequireITSStandAlone = 0;
333 fCutRequireITSpureSA = 0;
335 fCutNsigmaToVertex = 0;
336 fCutSigmaToVertexRequired = 0;
337 fCutMaxDCAToVertexXY = 0;
338 fCutMaxDCAToVertexZ = 0;
339 fCutDCAToVertex2D = 0;
340 fCutMinDCAToVertexXY = 0;
341 fCutMinDCAToVertexZ = 0;
342 fCutMaxDCAToVertexXYPtDep = "";
343 fCutMaxDCAToVertexZPtDep = "";
344 fCutMinDCAToVertexXYPtDep = "";
345 fCutMinDCAToVertexZPtDep = "";
347 if(f1CutMaxDCAToVertexXYPtDep)delete f1CutMaxDCAToVertexXYPtDep;
348 f1CutMaxDCAToVertexXYPtDep = 0;
349 if( f1CutMaxDCAToVertexXYPtDep) delete f1CutMaxDCAToVertexXYPtDep;
350 f1CutMaxDCAToVertexXYPtDep = 0;
351 if( f1CutMaxDCAToVertexZPtDep) delete f1CutMaxDCAToVertexZPtDep;
352 f1CutMaxDCAToVertexZPtDep = 0;
353 if( f1CutMinDCAToVertexXYPtDep)delete f1CutMinDCAToVertexXYPtDep;
354 f1CutMinDCAToVertexXYPtDep = 0;
355 if(f1CutMinDCAToVertexZPtDep)delete f1CutMinDCAToVertexZPtDep;
356 f1CutMinDCAToVertexZPtDep = 0;
374 fHistogramsOn = kFALSE;
376 for (Int_t i=0; i<2; ++i)
378 fhNClustersITS[i] = 0;
379 fhNClustersTPC[i] = 0;
381 fhChi2PerClusterITS[i] = 0;
382 fhChi2PerClusterTPC[i] = 0;
390 fhRel1PtUncertainty[i] = 0;
397 fhDXYNormalized[i] = 0;
398 fhDZNormalized[i] = 0;
399 fhDXYvsDZNormalized[i] = 0;
400 fhNSigmaToVertex[i] = 0;
408 fhCutCorrelation = 0;
411 //_____________________________________________________________________________
412 AliESDtrackCuts &AliESDtrackCuts::operator=(const AliESDtrackCuts &c)
415 // Assignment operator
418 if (this != &c) ((AliESDtrackCuts &) c).Copy(*this);
422 //_____________________________________________________________________________
423 void AliESDtrackCuts::Copy(TObject &c) const
429 AliESDtrackCuts& target = (AliESDtrackCuts &) c;
433 target.fCutMinNClusterTPC = fCutMinNClusterTPC;
434 target.fCutMinNClusterITS = fCutMinNClusterITS;
436 target.fCutMaxChi2PerClusterTPC = fCutMaxChi2PerClusterTPC;
437 target.fCutMaxChi2PerClusterITS = fCutMaxChi2PerClusterITS;
439 for (Int_t i = 0; i < 3; i++)
440 target.fCutClusterRequirementITS[i] = fCutClusterRequirementITS[i];
442 target.fCutMaxC11 = fCutMaxC11;
443 target.fCutMaxC22 = fCutMaxC22;
444 target.fCutMaxC33 = fCutMaxC33;
445 target.fCutMaxC44 = fCutMaxC44;
446 target.fCutMaxC55 = fCutMaxC55;
448 target.fCutMaxRel1PtUncertainty = fCutMaxRel1PtUncertainty;
450 target.fCutAcceptKinkDaughters = fCutAcceptKinkDaughters;
451 target.fCutAcceptSharedTPCClusters = fCutAcceptSharedTPCClusters;
452 target.fCutMaxFractionSharedTPCClusters = fCutMaxFractionSharedTPCClusters;
453 target.fCutRequireTPCRefit = fCutRequireTPCRefit;
454 target.fCutRequireTPCStandAlone = fCutRequireTPCStandAlone;
455 target.fCutRequireITSRefit = fCutRequireITSRefit;
456 target.fCutRequireITSPid = fCutRequireITSPid;
457 target.fCutRequireITSStandAlone = fCutRequireITSStandAlone;
458 target.fCutRequireITSpureSA = fCutRequireITSpureSA;
460 target.fCutNsigmaToVertex = fCutNsigmaToVertex;
461 target.fCutSigmaToVertexRequired = fCutSigmaToVertexRequired;
462 target.fCutMaxDCAToVertexXY = fCutMaxDCAToVertexXY;
463 target.fCutMaxDCAToVertexZ = fCutMaxDCAToVertexZ;
464 target.fCutDCAToVertex2D = fCutDCAToVertex2D;
465 target.fCutMinDCAToVertexXY = fCutMinDCAToVertexXY;
466 target.fCutMinDCAToVertexZ = fCutMinDCAToVertexZ;
468 target.fCutMaxDCAToVertexXYPtDep = fCutMaxDCAToVertexXYPtDep;
469 target.SetMaxDCAToVertexXYPtDep(fCutMaxDCAToVertexXYPtDep.Data());
471 target.fCutMaxDCAToVertexZPtDep = fCutMaxDCAToVertexZPtDep;
472 target.SetMaxDCAToVertexZPtDep(fCutMaxDCAToVertexZPtDep.Data());
474 target.fCutMinDCAToVertexXYPtDep = fCutMinDCAToVertexXYPtDep;
475 target.SetMinDCAToVertexXYPtDep(fCutMinDCAToVertexXYPtDep.Data());
477 target.fCutMinDCAToVertexZPtDep = fCutMinDCAToVertexZPtDep;
478 target.SetMinDCAToVertexZPtDep(fCutMinDCAToVertexZPtDep.Data());
480 target.fPMin = fPMin;
481 target.fPMax = fPMax;
482 target.fPtMin = fPtMin;
483 target.fPtMax = fPtMax;
484 target.fPxMin = fPxMin;
485 target.fPxMax = fPxMax;
486 target.fPyMin = fPyMin;
487 target.fPyMax = fPyMax;
488 target.fPzMin = fPzMin;
489 target.fPzMax = fPzMax;
490 target.fEtaMin = fEtaMin;
491 target.fEtaMax = fEtaMax;
492 target.fRapMin = fRapMin;
493 target.fRapMax = fRapMax;
495 target.fHistogramsOn = fHistogramsOn;
497 for (Int_t i=0; i<2; ++i)
499 if (fhNClustersITS[i]) target.fhNClustersITS[i] = (TH1F*) fhNClustersITS[i]->Clone();
500 if (fhNClustersTPC[i]) target.fhNClustersTPC[i] = (TH1F*) fhNClustersTPC[i]->Clone();
502 if (fhChi2PerClusterITS[i]) target.fhChi2PerClusterITS[i] = (TH1F*) fhChi2PerClusterITS[i]->Clone();
503 if (fhChi2PerClusterTPC[i]) target.fhChi2PerClusterTPC[i] = (TH1F*) fhChi2PerClusterTPC[i]->Clone();
505 if (fhC11[i]) target.fhC11[i] = (TH1F*) fhC11[i]->Clone();
506 if (fhC22[i]) target.fhC22[i] = (TH1F*) fhC22[i]->Clone();
507 if (fhC33[i]) target.fhC33[i] = (TH1F*) fhC33[i]->Clone();
508 if (fhC44[i]) target.fhC44[i] = (TH1F*) fhC44[i]->Clone();
509 if (fhC55[i]) target.fhC55[i] = (TH1F*) fhC55[i]->Clone();
511 if (fhRel1PtUncertainty[i]) target.fhRel1PtUncertainty[i] = (TH1F*) fhRel1PtUncertainty[i]->Clone();
513 if (fhDXY[i]) target.fhDXY[i] = (TH1F*) fhDXY[i]->Clone();
514 if (fhDZ[i]) target.fhDZ[i] = (TH1F*) fhDZ[i]->Clone();
515 if (fhDXYDZ[i]) target.fhDXYDZ[i] = (TH1F*) fhDXYDZ[i]->Clone();
516 if (fhDXYvsDZ[i]) target.fhDXYvsDZ[i] = (TH2F*) fhDXYvsDZ[i]->Clone();
518 if (fhDXYNormalized[i]) target.fhDXYNormalized[i] = (TH1F*) fhDXYNormalized[i]->Clone();
519 if (fhDZNormalized[i]) target.fhDZNormalized[i] = (TH1F*) fhDZNormalized[i]->Clone();
520 if (fhDXYvsDZNormalized[i]) target.fhDXYvsDZNormalized[i] = (TH2F*) fhDXYvsDZNormalized[i]->Clone();
521 if (fhNSigmaToVertex[i]) target.fhNSigmaToVertex[i] = (TH1F*) fhNSigmaToVertex[i]->Clone();
523 if (fhPt[i]) target.fhPt[i] = (TH1F*) fhPt[i]->Clone();
524 if (fhEta[i]) target.fhEta[i] = (TH1F*) fhEta[i]->Clone();
526 if (ffDTheoretical) target.ffDTheoretical = (TF1*) ffDTheoretical->Clone();
528 if (fhCutStatistics) target.fhCutStatistics = (TH1F*) fhCutStatistics->Clone();
529 if (fhCutCorrelation) target.fhCutCorrelation = (TH2F*) fhCutCorrelation->Clone();
534 //_____________________________________________________________________________
535 Long64_t AliESDtrackCuts::Merge(TCollection* list) {
536 // Merge a list of AliESDtrackCuts objects with this (needed for PROOF)
537 // Returns the number of merged objects (including this)
544 TIterator* iter = list->MakeIterator();
547 // collection of measured and generated histograms
549 while ((obj = iter->Next())) {
551 AliESDtrackCuts* entry = dynamic_cast<AliESDtrackCuts*>(obj);
555 if (!entry->fHistogramsOn)
558 for (Int_t i=0; i<2; i++) {
560 fhNClustersITS[i] ->Add(entry->fhNClustersITS[i] );
561 fhNClustersTPC[i] ->Add(entry->fhNClustersTPC[i] );
563 fhChi2PerClusterITS[i] ->Add(entry->fhChi2PerClusterITS[i]);
564 fhChi2PerClusterTPC[i] ->Add(entry->fhChi2PerClusterTPC[i]);
566 fhC11[i] ->Add(entry->fhC11[i] );
567 fhC22[i] ->Add(entry->fhC22[i] );
568 fhC33[i] ->Add(entry->fhC33[i] );
569 fhC44[i] ->Add(entry->fhC44[i] );
570 fhC55[i] ->Add(entry->fhC55[i] );
572 fhRel1PtUncertainty[i] ->Add(entry->fhRel1PtUncertainty[i]);
574 fhDXY[i] ->Add(entry->fhDXY[i] );
575 fhDZ[i] ->Add(entry->fhDZ[i] );
576 fhDXYDZ[i] ->Add(entry->fhDXYDZ[i] );
577 fhDXYvsDZ[i] ->Add(entry->fhDXYvsDZ[i] );
579 fhDXYNormalized[i] ->Add(entry->fhDXYNormalized[i] );
580 fhDZNormalized[i] ->Add(entry->fhDZNormalized[i] );
581 fhDXYvsDZNormalized[i] ->Add(entry->fhDXYvsDZNormalized[i]);
582 fhNSigmaToVertex[i] ->Add(entry->fhNSigmaToVertex[i]);
584 fhPt[i] ->Add(entry->fhPt[i]);
585 fhEta[i] ->Add(entry->fhEta[i]);
588 fhCutStatistics ->Add(entry->fhCutStatistics);
589 fhCutCorrelation ->Add(entry->fhCutCorrelation);
596 //____________________________________________________________________
597 AliESDtrackCuts* AliESDtrackCuts::GetStandardTPCOnlyTrackCuts()
599 // creates an AliESDtrackCuts object and fills it with standard (pre data-taking) values for TPC-only cuts
601 Printf("AliESDtrackCuts::GetStandardTPCOnlyTrackCuts: Creating track cuts for TPC-only.");
603 AliESDtrackCuts* esdTrackCuts = new AliESDtrackCuts;
605 esdTrackCuts->SetMinNClustersTPC(50);
606 esdTrackCuts->SetMaxChi2PerClusterTPC(4);
607 esdTrackCuts->SetAcceptKinkDaughters(kFALSE);
609 esdTrackCuts->SetMaxDCAToVertexZ(3.2);
610 esdTrackCuts->SetMaxDCAToVertexXY(2.4);
611 esdTrackCuts->SetDCAToVertex2D(kTRUE);
616 //____________________________________________________________________
617 AliESDtrackCuts* AliESDtrackCuts::GetStandardITSTPCTrackCuts2009(Bool_t selPrimaries)
619 // creates an AliESDtrackCuts object and fills it with standard values for ITS-TPC cuts for pp 2009 data
621 Printf("AliESDtrackCuts::GetStandardITSTPCTrackCuts: Creating track cuts for ITS+TPC.");
623 AliESDtrackCuts* esdTrackCuts = new AliESDtrackCuts;
626 esdTrackCuts->SetRequireTPCStandAlone(kTRUE); // to get chi2 and ncls of kTPCin
627 esdTrackCuts->SetMinNClustersTPC(70);
628 esdTrackCuts->SetMaxChi2PerClusterTPC(4);
629 esdTrackCuts->SetAcceptKinkDaughters(kFALSE);
630 esdTrackCuts->SetRequireTPCRefit(kTRUE);
632 esdTrackCuts->SetRequireITSRefit(kTRUE);
633 esdTrackCuts->SetClusterRequirementITS(AliESDtrackCuts::kSPD,
634 AliESDtrackCuts::kAny);
636 // 7*(0.0050+0.0060/pt^0.9)
637 esdTrackCuts->SetMaxDCAToVertexXYPtDep("0.0350+0.0420/pt^0.9");
639 esdTrackCuts->SetMaxDCAToVertexZ(1.e6);
640 esdTrackCuts->SetDCAToVertex2D(kFALSE);
641 esdTrackCuts->SetRequireSigmaToVertex(kFALSE);
642 //esdTrackCuts->SetEtaRange(-0.8,+0.8);
647 //____________________________________________________________________
648 AliESDtrackCuts* AliESDtrackCuts::GetStandardITSTPCTrackCuts2010(Bool_t selPrimaries)
650 // creates an AliESDtrackCuts object and fills it with standard values for ITS-TPC cuts for pp 2010 data
652 Printf("AliESDtrackCuts::GetStandardITSTPCTrackCuts: Creating track cuts for ITS+TPC.");
654 AliESDtrackCuts* esdTrackCuts = new AliESDtrackCuts;
657 esdTrackCuts->SetMinNClustersTPC(70);
658 esdTrackCuts->SetMaxChi2PerClusterTPC(4);
659 esdTrackCuts->SetAcceptKinkDaughters(kFALSE);
660 esdTrackCuts->SetRequireTPCRefit(kTRUE);
662 esdTrackCuts->SetRequireITSRefit(kTRUE);
663 esdTrackCuts->SetClusterRequirementITS(AliESDtrackCuts::kSPD,
664 AliESDtrackCuts::kAny);
666 // 7*(0.0026+0.0050/pt^1.01)
667 esdTrackCuts->SetMaxDCAToVertexXYPtDep("0.0182+0.0350/pt^1.01");
669 esdTrackCuts->SetMaxDCAToVertexZ(2);
670 esdTrackCuts->SetDCAToVertex2D(kFALSE);
671 esdTrackCuts->SetRequireSigmaToVertex(kFALSE);
676 //____________________________________________________________________
677 AliESDtrackCuts* AliESDtrackCuts::GetStandardITSPureSATrackCuts2009(Bool_t selPrimaries, Bool_t useForPid)
679 // creates an AliESDtrackCuts object and fills it with standard values for ITS pure SA tracks
681 AliESDtrackCuts* esdTrackCuts = new AliESDtrackCuts;
682 esdTrackCuts->SetRequireITSPureStandAlone(kTRUE);
683 esdTrackCuts->SetRequireITSRefit(kTRUE);
684 esdTrackCuts->SetMinNClustersITS(4);
685 esdTrackCuts->SetClusterRequirementITS(AliESDtrackCuts::kSPD,
686 AliESDtrackCuts::kAny);
687 esdTrackCuts->SetMaxChi2PerClusterITS(1.);
690 // 7*(0.0085+0.0026/pt^1.55)
691 esdTrackCuts->SetMaxDCAToVertexXYPtDep("0.0595+0.0182/pt^1.55");
694 esdTrackCuts->SetRequireITSPid(kTRUE);
699 //____________________________________________________________________
700 AliESDtrackCuts* AliESDtrackCuts::GetStandardITSSATrackCuts2009(Bool_t selPrimaries, Bool_t useForPid)
702 // creates an AliESDtrackCuts object and fills it with standard values for ITS pure SA tracks
704 AliESDtrackCuts* esdTrackCuts = new AliESDtrackCuts;
705 esdTrackCuts->SetRequireITSStandAlone(kTRUE);
706 esdTrackCuts->SetRequireITSPureStandAlone(kFALSE);
707 esdTrackCuts->SetRequireITSRefit(kTRUE);
708 esdTrackCuts->SetMinNClustersITS(4);
709 esdTrackCuts->SetClusterRequirementITS(AliESDtrackCuts::kSPD,
710 AliESDtrackCuts::kAny);
711 esdTrackCuts->SetMaxChi2PerClusterITS(1.);
714 // 7*(0.0085+0.0026/pt^1.55)
715 esdTrackCuts->SetMaxDCAToVertexXYPtDep("0.0595+0.0182/pt^1.55");
718 esdTrackCuts->SetRequireITSPid(kTRUE);
724 //____________________________________________________________________
725 Int_t AliESDtrackCuts::GetReferenceMultiplicity(AliESDEvent* esd, Bool_t tpcOnly)
727 // Gets reference multiplicity following the standard cuts and a defined fiducial volume
728 // tpcOnly = kTRUE -> consider TPC-only tracks
729 // = kFALSE -> consider global tracks
733 Printf("AliESDtrackCuts::GetReferenceMultiplicity: Not implemented for global tracks!");
737 AliESDtrackCuts* esdTrackCuts = GetStandardTPCOnlyTrackCuts();
738 esdTrackCuts->SetEtaRange(-0.8, 0.8);
739 esdTrackCuts->SetPtRange(0.15);
741 Int_t nTracks = esdTrackCuts->CountAcceptedTracks(esd);
749 //____________________________________________________________________
750 Float_t AliESDtrackCuts::GetSigmaToVertex(AliESDtrack* esdTrack)
752 // Calculates the number of sigma to the vertex.
757 esdTrack->GetImpactParameters(b,bCov);
759 if (bCov[0]<=0 || bCov[2]<=0) {
760 AliDebugClass(1, "Estimated b resolution lower or equal zero!");
761 bCov[0]=0; bCov[2]=0;
763 bRes[0] = TMath::Sqrt(bCov[0]);
764 bRes[1] = TMath::Sqrt(bCov[2]);
766 // -----------------------------------
767 // How to get to a n-sigma cut?
769 // The accumulated statistics from 0 to d is
771 // -> Erf(d/Sqrt(2)) for a 1-dim gauss (d = n_sigma)
772 // -> 1 - Exp(-d**2) for a 2-dim gauss (d*d = dx*dx + dy*dy != n_sigma)
774 // It means that for a 2-dim gauss: n_sigma(d) = Sqrt(2)*ErfInv(1 - Exp((-d**2)/2)
775 // Can this be expressed in a different way?
777 if (bRes[0] == 0 || bRes[1] ==0)
780 Float_t d = TMath::Sqrt(TMath::Power(b[0]/bRes[0],2) + TMath::Power(b[1]/bRes[1],2));
782 // work around precision problem
783 // if d is too big, TMath::Exp(...) gets 0, and TMath::ErfInverse(1) that should be infinite, gets 0 :(
784 // 1e-15 corresponds to nsigma ~ 7.7
785 if (TMath::Exp(-d * d / 2) < 1e-15)
788 Float_t nSigma = TMath::ErfInverse(1 - TMath::Exp(-d * d / 2)) * TMath::Sqrt(2);
792 void AliESDtrackCuts::EnableNeededBranches(TTree* tree)
794 // enables the branches needed by AcceptTrack, for a list see comment of AcceptTrack
796 tree->SetBranchStatus("fTracks.fFlags", 1);
797 tree->SetBranchStatus("fTracks.fITSncls", 1);
798 tree->SetBranchStatus("fTracks.fTPCncls", 1);
799 tree->SetBranchStatus("fTracks.fITSchi2", 1);
800 tree->SetBranchStatus("fTracks.fTPCchi2", 1);
801 tree->SetBranchStatus("fTracks.fC*", 1);
802 tree->SetBranchStatus("fTracks.fD", 1);
803 tree->SetBranchStatus("fTracks.fZ", 1);
804 tree->SetBranchStatus("fTracks.fCdd", 1);
805 tree->SetBranchStatus("fTracks.fCdz", 1);
806 tree->SetBranchStatus("fTracks.fCzz", 1);
807 tree->SetBranchStatus("fTracks.fP*", 1);
808 tree->SetBranchStatus("fTracks.fR*", 1);
809 tree->SetBranchStatus("fTracks.fKinkIndexes*", 1);
812 //____________________________________________________________________
813 Bool_t AliESDtrackCuts::AcceptTrack(AliESDtrack* esdTrack)
816 // figure out if the tracks survives all the track cuts defined
818 // the different quality parameter and kinematic values are first
819 // retrieved from the track. then it is found out what cuts the
820 // track did not survive and finally the cuts are imposed.
822 // this function needs the following branches:
828 // fTracks.fC //GetExternalCovariance
829 // fTracks.fD //GetImpactParameters
830 // fTracks.fZ //GetImpactParameters
831 // fTracks.fCdd //GetImpactParameters
832 // fTracks.fCdz //GetImpactParameters
833 // fTracks.fCzz //GetImpactParameters
834 // fTracks.fP //GetPxPyPz
835 // fTracks.fR //GetMass
836 // fTracks.fP //GetMass
837 // fTracks.fKinkIndexes
839 UInt_t status = esdTrack->GetStatus();
841 // getting quality parameters from the ESD track
842 Int_t nClustersITS = esdTrack->GetITSclusters(0);
843 Int_t nClustersTPC = -1;
844 if(fCutRequireTPCStandAlone) {
845 nClustersTPC = esdTrack->GetTPCNclsIter1();
848 nClustersTPC = esdTrack->GetTPCclusters(0);
851 Int_t nClustersTPCShared = esdTrack->GetTPCnclsS();
852 Float_t fracClustersTPCShared = -1.;
854 Float_t chi2PerClusterITS = -1;
855 Float_t chi2PerClusterTPC = -1;
857 chi2PerClusterITS = esdTrack->GetITSchi2()/Float_t(nClustersITS);
858 if (nClustersTPC!=0) {
859 if(fCutRequireTPCStandAlone) {
860 chi2PerClusterTPC = esdTrack->GetTPCchi2Iter1()/Float_t(nClustersTPC);
862 chi2PerClusterTPC = esdTrack->GetTPCchi2()/Float_t(nClustersTPC);
864 fracClustersTPCShared = Float_t(nClustersTPCShared)/Float_t(nClustersTPC);
868 esdTrack->GetExternalCovariance(extCov);
870 // getting the track to vertex parameters
871 Float_t nSigmaToVertex = GetSigmaToVertex(esdTrack);
875 esdTrack->GetImpactParameters(b,bCov);
876 if (bCov[0]<=0 || bCov[2]<=0) {
877 AliDebug(1, "Estimated b resolution lower or equal zero!");
878 bCov[0]=0; bCov[2]=0;
882 // set pt-dependent DCA cuts, if requested
883 SetPtDepDCACuts(esdTrack->Pt());
886 Float_t dcaToVertexXY = b[0];
887 Float_t dcaToVertexZ = b[1];
889 Float_t dcaToVertex = -1;
891 if (fCutDCAToVertex2D)
893 dcaToVertex = TMath::Sqrt(dcaToVertexXY*dcaToVertexXY/fCutMaxDCAToVertexXY/fCutMaxDCAToVertexXY + dcaToVertexZ*dcaToVertexZ/fCutMaxDCAToVertexZ/fCutMaxDCAToVertexZ);
896 dcaToVertex = TMath::Sqrt(dcaToVertexXY*dcaToVertexXY + dcaToVertexZ*dcaToVertexZ);
898 // getting the kinematic variables of the track
899 // (assuming the mass is known)
901 esdTrack->GetPxPyPz(p);
903 Float_t momentum = TMath::Sqrt(TMath::Power(p[0],2) + TMath::Power(p[1],2) + TMath::Power(p[2],2));
904 Float_t pt = TMath::Sqrt(TMath::Power(p[0],2) + TMath::Power(p[1],2));
905 Float_t energy = TMath::Sqrt(TMath::Power(esdTrack->GetMass(),2) + TMath::Power(momentum,2));
907 //y-eta related calculations
910 if((momentum != TMath::Abs(p[2]))&&(momentum != 0))
911 eta = 0.5*TMath::Log((momentum + p[2])/(momentum - p[2]));
912 if((energy != TMath::Abs(p[2]))&&(momentum != 0))
913 y = 0.5*TMath::Log((energy + p[2])/(energy - p[2]));
917 Printf("AliESDtrackCuts::AcceptTrack: WARNING: GetSigma1Pt2() returns negative value for external covariance matrix element fC[14]: %f. Corrupted track information, track will not be accepted!", extCov[14]);
920 Float_t relUncertainty1Pt = TMath::Sqrt(extCov[14])*pt;
922 //########################################################################
926 for (Int_t i=0; i<kNCuts; i++) cuts[i]=kFALSE;
928 // track quality cuts
929 if (fCutRequireTPCRefit && (status&AliESDtrack::kTPCrefit)==0)
931 if (fCutRequireTPCStandAlone && (status&AliESDtrack::kTPCin)==0)
933 if (fCutRequireITSRefit && (status&AliESDtrack::kITSrefit)==0)
935 if (nClustersTPC<fCutMinNClusterTPC)
937 if (nClustersITS<fCutMinNClusterITS)
939 if (chi2PerClusterTPC>fCutMaxChi2PerClusterTPC)
941 if (chi2PerClusterITS>fCutMaxChi2PerClusterITS)
943 if (extCov[0] > fCutMaxC11)
945 if (extCov[2] > fCutMaxC22)
947 if (extCov[5] > fCutMaxC33)
949 if (extCov[9] > fCutMaxC44)
951 if (extCov[14] > fCutMaxC55)
953 if (nSigmaToVertex > fCutNsigmaToVertex && fCutSigmaToVertexRequired)
955 // if n sigma could not be calculated
956 if (nSigmaToVertex<0 && fCutSigmaToVertexRequired)
958 if (!fCutAcceptKinkDaughters && esdTrack->GetKinkIndex(0)>0)
960 // track kinematics cut
961 if((momentum < fPMin) || (momentum > fPMax))
963 if((pt < fPtMin) || (pt > fPtMax))
965 if((p[0] < fPxMin) || (p[0] > fPxMax))
967 if((p[1] < fPyMin) || (p[1] > fPyMax))
969 if((p[2] < fPzMin) || (p[2] > fPzMax))
971 if((eta < fEtaMin) || (eta > fEtaMax))
973 if((y < fRapMin) || (y > fRapMax))
975 if (fCutDCAToVertex2D && dcaToVertex > 1)
977 if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexXY) > fCutMaxDCAToVertexXY)
979 if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexZ) > fCutMaxDCAToVertexZ)
981 if (fCutDCAToVertex2D && fCutMinDCAToVertexXY > 0 && fCutMinDCAToVertexZ > 0 && dcaToVertexXY*dcaToVertexXY/fCutMinDCAToVertexXY/fCutMinDCAToVertexXY + dcaToVertexZ*dcaToVertexZ/fCutMinDCAToVertexZ/fCutMinDCAToVertexZ < 1)
983 if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexXY) < fCutMinDCAToVertexXY)
985 if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexZ) < fCutMinDCAToVertexZ)
988 for (Int_t i = 0; i < 3; i++)
989 cuts[28+i] = !CheckITSClusterRequirement(fCutClusterRequirementITS[i], esdTrack->HasPointOnITSLayer(i*2), esdTrack->HasPointOnITSLayer(i*2+1));
991 if(fCutRequireITSStandAlone || fCutRequireITSpureSA){
992 if ((status & AliESDtrack::kITSin) == 0 || (status & AliESDtrack::kTPCin)){
996 // ITS standalone tracks
997 if(fCutRequireITSStandAlone && !fCutRequireITSpureSA){
998 if(status & AliESDtrack::kITSpureSA) cuts[31] = kTRUE;
999 }else if(fCutRequireITSpureSA){
1000 if(!(status & AliESDtrack::kITSpureSA)) cuts[31] = kTRUE;
1005 if (relUncertainty1Pt > fCutMaxRel1PtUncertainty)
1008 if (!fCutAcceptSharedTPCClusters && nClustersTPCShared!=0)
1011 if (fracClustersTPCShared > fCutMaxFractionSharedTPCClusters)
1014 if(fCutRequireITSPid){
1015 UChar_t clumap=esdTrack->GetITSClusterMap();
1016 Int_t nPointsForPid=0;
1017 for(Int_t i=2; i<6; i++){
1018 if(clumap&(1<<i)) ++nPointsForPid;
1020 if(nPointsForPid<3) cuts[35] = kTRUE;
1024 for (Int_t i=0; i<kNCuts; i++)
1025 if (cuts[i]) {cut = kTRUE;}
1028 //########################################################################
1029 // filling histograms
1030 if (fHistogramsOn) {
1031 fhCutStatistics->Fill(fhCutStatistics->GetBinCenter(fhCutStatistics->GetXaxis()->FindBin("n tracks")));
1033 fhCutStatistics->Fill(fhCutStatistics->GetBinCenter(fhCutStatistics->GetXaxis()->FindBin("n cut tracks")));
1035 for (Int_t i=0; i<kNCuts; i++) {
1037 fhCutStatistics->Fill(fhCutStatistics->GetBinCenter(fhCutStatistics->GetXaxis()->FindBin(fgkCutNames[i])));
1039 for (Int_t j=i; j<kNCuts; j++) {
1040 if (cuts[i] && cuts[j]) {
1041 Float_t xC = fhCutCorrelation->GetXaxis()->GetBinCenter(fhCutCorrelation->GetXaxis()->FindBin(fgkCutNames[i]));
1042 Float_t yC = fhCutCorrelation->GetYaxis()->GetBinCenter(fhCutCorrelation->GetYaxis()->FindBin(fgkCutNames[j]));
1043 fhCutCorrelation->Fill(xC, yC);
1049 // now we loop over the filling of the histograms twice: once "before" the cut, once "after"
1050 // the code is not in a function due to too many local variables that would need to be passed
1052 for (Int_t id = 0; id < 2; id++)
1054 // id = 0 --> before cut
1055 // id = 1 --> after cut
1059 fhNClustersITS[id]->Fill(nClustersITS);
1060 fhNClustersTPC[id]->Fill(nClustersTPC);
1061 fhChi2PerClusterITS[id]->Fill(chi2PerClusterITS);
1062 fhChi2PerClusterTPC[id]->Fill(chi2PerClusterTPC);
1064 fhC11[id]->Fill(extCov[0]);
1065 fhC22[id]->Fill(extCov[2]);
1066 fhC33[id]->Fill(extCov[5]);
1067 fhC44[id]->Fill(extCov[9]);
1068 fhC55[id]->Fill(extCov[14]);
1070 fhRel1PtUncertainty[id]->Fill(relUncertainty1Pt);
1073 fhEta[id]->Fill(eta);
1076 bRes[0] = TMath::Sqrt(bCov[0]);
1077 bRes[1] = TMath::Sqrt(bCov[2]);
1079 fhDZ[id]->Fill(b[1]);
1080 fhDXY[id]->Fill(b[0]);
1081 fhDXYDZ[id]->Fill(dcaToVertex);
1082 fhDXYvsDZ[id]->Fill(b[1],b[0]);
1084 if (bRes[0]!=0 && bRes[1]!=0) {
1085 fhDZNormalized[id]->Fill(b[1]/bRes[1]);
1086 fhDXYNormalized[id]->Fill(b[0]/bRes[0]);
1087 fhDXYvsDZNormalized[id]->Fill(b[1]/bRes[1], b[0]/bRes[0]);
1088 fhNSigmaToVertex[id]->Fill(nSigmaToVertex);
1100 //____________________________________________________________________
1101 Bool_t AliESDtrackCuts::CheckITSClusterRequirement(ITSClusterRequirement req, Bool_t clusterL1, Bool_t clusterL2)
1103 // checks if the cluster requirement is fullfilled (in this case: return kTRUE)
1107 case kOff: return kTRUE;
1108 case kNone: return !clusterL1 && !clusterL2;
1109 case kAny: return clusterL1 || clusterL2;
1110 case kFirst: return clusterL1;
1111 case kOnlyFirst: return clusterL1 && !clusterL2;
1112 case kSecond: return clusterL2;
1113 case kOnlySecond: return clusterL2 && !clusterL1;
1114 case kBoth: return clusterL1 && clusterL2;
1120 //____________________________________________________________________
1121 AliESDtrack* AliESDtrackCuts::GetTPCOnlyTrack(AliESDEvent* esd, Int_t iTrack)
1124 // Utility function to
1125 // create a TPC only track from the given esd track
1127 // IMPORTANT: The track has to be deleted by the user
1129 // NB. most of the functionality to get a TPC only track from an ESD track is in AliESDtrack, where it should be
1130 // there are only missing propagations here that are needed for old data
1131 // this function will therefore become obsolete
1133 // adapted from code provided by CKB
1135 if (!esd->GetPrimaryVertexTPC())
1136 return 0; // No TPC vertex no TPC tracks
1138 if(!esd->GetPrimaryVertexTPC()->GetStatus())
1139 return 0; // TPC Vertex is created by default in AliESDEvent, do not use in this case
1141 AliESDtrack* track = esd->GetTrack(iTrack);
1145 AliESDtrack *tpcTrack = new AliESDtrack();
1147 // only true if we have a tpc track
1148 if (!track->FillTPCOnlyTrack(*tpcTrack))
1154 // propagate to Vertex
1155 // not needed for normal reconstructed ESDs...
1156 // Double_t pTPC[2],covTPC[3];
1157 // tpcTrack->PropagateToDCA(esd->GetPrimaryVertexTPC(), esd->GetMagneticField(), 10000, pTPC, covTPC);
1162 //____________________________________________________________________
1163 TObjArray* AliESDtrackCuts::GetAcceptedTracks(AliESDEvent* esd,Bool_t bTPC)
1166 // returns an array of all tracks that pass the cuts
1167 // or an array of TPC only tracks (propagated to the TPC vertex during reco)
1168 // tracks that pass the cut
1170 TObjArray* acceptedTracks = new TObjArray();
1172 // loop over esd tracks
1173 for (Int_t iTrack = 0; iTrack < esd->GetNumberOfTracks(); iTrack++) {
1175 if(!esd->GetPrimaryVertexTPC())return acceptedTracks; // No TPC vertex no TPC tracks
1176 if(!esd->GetPrimaryVertexTPC()->GetStatus())return acceptedTracks; // No proper TPC vertex, only the default
1178 AliESDtrack *tpcTrack = GetTPCOnlyTrack(esd, iTrack);
1182 if (AcceptTrack(tpcTrack)) {
1183 acceptedTracks->Add(tpcTrack);
1190 AliESDtrack* track = esd->GetTrack(iTrack);
1191 if(AcceptTrack(track))
1192 acceptedTracks->Add(track);
1195 if(bTPC)acceptedTracks->SetOwner(kTRUE);
1196 return acceptedTracks;
1199 //____________________________________________________________________
1200 Int_t AliESDtrackCuts::CountAcceptedTracks(AliESDEvent* esd)
1203 // returns an the number of tracks that pass the cuts
1208 // loop over esd tracks
1209 for (Int_t iTrack = 0; iTrack < esd->GetNumberOfTracks(); iTrack++) {
1210 AliESDtrack* track = esd->GetTrack(iTrack);
1211 if (AcceptTrack(track))
1218 //____________________________________________________________________
1219 void AliESDtrackCuts::DefineHistograms(Int_t color) {
1221 // diagnostics histograms are defined
1224 fHistogramsOn=kTRUE;
1226 Bool_t oldStatus = TH1::AddDirectoryStatus();
1227 TH1::AddDirectory(kFALSE);
1229 //###################################################################################
1230 // defining histograms
1232 fhCutStatistics = new TH1F("cut_statistics","cut statistics",kNCuts+4,-0.5,kNCuts+3.5);
1234 fhCutStatistics->GetXaxis()->SetBinLabel(1,"n tracks");
1235 fhCutStatistics->GetXaxis()->SetBinLabel(2,"n cut tracks");
1237 fhCutCorrelation = new TH2F("cut_correlation","cut correlation",kNCuts,-0.5,kNCuts-0.5,kNCuts,-0.5,kNCuts-0.5);;
1239 for (Int_t i=0; i<kNCuts; i++) {
1240 fhCutStatistics->GetXaxis()->SetBinLabel(i+4,fgkCutNames[i]);
1241 fhCutCorrelation->GetXaxis()->SetBinLabel(i+1,fgkCutNames[i]);
1242 fhCutCorrelation->GetYaxis()->SetBinLabel(i+1,fgkCutNames[i]);
1245 fhCutStatistics ->SetLineColor(color);
1246 fhCutCorrelation ->SetLineColor(color);
1247 fhCutStatistics ->SetLineWidth(2);
1248 fhCutCorrelation ->SetLineWidth(2);
1250 for (Int_t i=0; i<2; i++) {
1251 fhNClustersITS[i] = new TH1F("nClustersITS" ,"",8,-0.5,7.5);
1252 fhNClustersTPC[i] = new TH1F("nClustersTPC" ,"",165,-0.5,164.5);
1253 fhChi2PerClusterITS[i] = new TH1F("chi2PerClusterITS","",500,0,10);
1254 fhChi2PerClusterTPC[i] = new TH1F("chi2PerClusterTPC","",500,0,10);
1256 fhC11[i] = new TH1F("covMatrixDiagonal11","",2000,0,20);
1257 fhC22[i] = new TH1F("covMatrixDiagonal22","",2000,0,20);
1258 fhC33[i] = new TH1F("covMatrixDiagonal33","",1000,0,0.1);
1259 fhC44[i] = new TH1F("covMatrixDiagonal44","",1000,0,0.1);
1260 fhC55[i] = new TH1F("covMatrixDiagonal55","",1000,0,5);
1262 fhRel1PtUncertainty[i] = new TH1F("rel1PtUncertainty","",1000,0,5);
1264 fhDXY[i] = new TH1F("dXY" ,"",500,-10,10);
1265 fhDZ[i] = new TH1F("dZ" ,"",500,-10,10);
1266 fhDXYDZ[i] = new TH1F("dXYDZ" ,"",500,0,10);
1267 fhDXYvsDZ[i] = new TH2F("dXYvsDZ","",200,-10,10,200,-10,10);
1269 fhDXYNormalized[i] = new TH1F("dXYNormalized" ,"",500,-10,10);
1270 fhDZNormalized[i] = new TH1F("dZNormalized" ,"",500,-10,10);
1271 fhDXYvsDZNormalized[i] = new TH2F("dXYvsDZNormalized","",200,-10,10,200,-10,10);
1273 fhNSigmaToVertex[i] = new TH1F("nSigmaToVertex","",500,0,10);
1275 fhPt[i] = new TH1F("pt" ,"p_{T} distribution;p_{T} (GeV/c)", 800, 0.0, 10.0);
1276 fhEta[i] = new TH1F("eta" ,"#eta distribution;#eta",40,-2.0,2.0);
1278 fhNClustersITS[i]->SetTitle("n ITS clusters");
1279 fhNClustersTPC[i]->SetTitle("n TPC clusters");
1280 fhChi2PerClusterITS[i]->SetTitle("#Chi^{2} per ITS cluster");
1281 fhChi2PerClusterTPC[i]->SetTitle("#Chi^{2} per TPC cluster");
1283 fhC11[i]->SetTitle("cov 11 : #sigma_{y}^{2} [cm^{2}]");
1284 fhC22[i]->SetTitle("cov 22 : #sigma_{z}^{2} [cm^{2}]");
1285 fhC33[i]->SetTitle("cov 33 : #sigma_{sin(#phi)}^{2}");
1286 fhC44[i]->SetTitle("cov 44 : #sigma_{tan(#theta_{dip})}^{2}");
1287 fhC55[i]->SetTitle("cov 55 : #sigma_{1/p_{T}}^{2} [(c/GeV)^2]");
1289 fhRel1PtUncertainty[i]->SetTitle("rel. uncertainty of 1/p_{T}");
1291 fhDXY[i]->SetXTitle("transverse impact parameter (cm)");
1292 fhDZ[i]->SetXTitle("longitudinal impact parameter (cm)");
1293 fhDXYDZ[i]->SetTitle("absolute impact parameter;sqrt(dXY**2 + dZ**2) (cm)");
1294 fhDXYvsDZ[i]->SetXTitle("longitudinal impact parameter (cm)");
1295 fhDXYvsDZ[i]->SetYTitle("transverse impact parameter (cm)");
1297 fhDXYNormalized[i]->SetTitle("normalized trans impact par (n#sigma)");
1298 fhDZNormalized[i]->SetTitle("normalized long impact par (n#sigma)");
1299 fhDXYvsDZNormalized[i]->SetTitle("normalized long impact par (n#sigma)");
1300 fhDXYvsDZNormalized[i]->SetYTitle("normalized trans impact par (n#sigma)");
1301 fhNSigmaToVertex[i]->SetTitle("n #sigma to vertex");
1303 fhNClustersITS[i]->SetLineColor(color); fhNClustersITS[i]->SetLineWidth(2);
1304 fhNClustersTPC[i]->SetLineColor(color); fhNClustersTPC[i]->SetLineWidth(2);
1305 fhChi2PerClusterITS[i]->SetLineColor(color); fhChi2PerClusterITS[i]->SetLineWidth(2);
1306 fhChi2PerClusterTPC[i]->SetLineColor(color); fhChi2PerClusterTPC[i]->SetLineWidth(2);
1308 fhC11[i]->SetLineColor(color); fhC11[i]->SetLineWidth(2);
1309 fhC22[i]->SetLineColor(color); fhC22[i]->SetLineWidth(2);
1310 fhC33[i]->SetLineColor(color); fhC33[i]->SetLineWidth(2);
1311 fhC44[i]->SetLineColor(color); fhC44[i]->SetLineWidth(2);
1312 fhC55[i]->SetLineColor(color); fhC55[i]->SetLineWidth(2);
1314 fhRel1PtUncertainty[i]->SetLineColor(color); fhRel1PtUncertainty[i]->SetLineWidth(2);
1316 fhDXY[i]->SetLineColor(color); fhDXY[i]->SetLineWidth(2);
1317 fhDZ[i]->SetLineColor(color); fhDZ[i]->SetLineWidth(2);
1318 fhDXYDZ[i]->SetLineColor(color); fhDXYDZ[i]->SetLineWidth(2);
1320 fhDXYNormalized[i]->SetLineColor(color); fhDXYNormalized[i]->SetLineWidth(2);
1321 fhDZNormalized[i]->SetLineColor(color); fhDZNormalized[i]->SetLineWidth(2);
1322 fhNSigmaToVertex[i]->SetLineColor(color); fhNSigmaToVertex[i]->SetLineWidth(2);
1325 // The number of sigmas to the vertex is per definition gaussian
1326 ffDTheoretical = new TF1("nSigmaToVertexTheoretical","([0]/2.506628274)*exp(-(x**2)/2)",0,50);
1327 ffDTheoretical->SetParameter(0,1);
1329 TH1::AddDirectory(oldStatus);
1332 //____________________________________________________________________
1333 Bool_t AliESDtrackCuts::LoadHistograms(const Char_t* dir)
1336 // loads the histograms from a file
1337 // if dir is empty a directory with the name of this object is taken (like in SaveHistogram)
1343 if (!gDirectory->cd(dir))
1346 ffDTheoretical = dynamic_cast<TF1*> (gDirectory->Get("nSigmaToVertexTheory"));
1348 fhCutStatistics = dynamic_cast<TH1F*> (gDirectory->Get("cut_statistics"));
1349 fhCutCorrelation = dynamic_cast<TH2F*> (gDirectory->Get("cut_correlation"));
1351 for (Int_t i=0; i<2; i++) {
1354 gDirectory->cd("before_cuts");
1357 gDirectory->cd("after_cuts");
1359 fhNClustersITS[i] = dynamic_cast<TH1F*> (gDirectory->Get("nClustersITS" ));
1360 fhNClustersTPC[i] = dynamic_cast<TH1F*> (gDirectory->Get("nClustersTPC" ));
1361 fhChi2PerClusterITS[i] = dynamic_cast<TH1F*> (gDirectory->Get("chi2PerClusterITS"));
1362 fhChi2PerClusterTPC[i] = dynamic_cast<TH1F*> (gDirectory->Get("chi2PerClusterTPC"));
1364 fhC11[i] = dynamic_cast<TH1F*> (gDirectory->Get("covMatrixDiagonal11"));
1365 fhC22[i] = dynamic_cast<TH1F*> (gDirectory->Get("covMatrixDiagonal22"));
1366 fhC33[i] = dynamic_cast<TH1F*> (gDirectory->Get("covMatrixDiagonal33"));
1367 fhC44[i] = dynamic_cast<TH1F*> (gDirectory->Get("covMatrixDiagonal44"));
1368 fhC55[i] = dynamic_cast<TH1F*> (gDirectory->Get("covMatrixDiagonal55"));
1370 fhRel1PtUncertainty[i] = dynamic_cast<TH1F*> (gDirectory->Get("rel1PtUncertainty"));
1372 fhDXY[i] = dynamic_cast<TH1F*> (gDirectory->Get("dXY" ));
1373 fhDZ[i] = dynamic_cast<TH1F*> (gDirectory->Get("dZ" ));
1374 fhDXYDZ[i] = dynamic_cast<TH1F*> (gDirectory->Get("dXYDZ"));
1375 fhDXYvsDZ[i] = dynamic_cast<TH2F*> (gDirectory->Get("dXYvsDZ"));
1377 fhDXYNormalized[i] = dynamic_cast<TH1F*> (gDirectory->Get("dXYNormalized" ));
1378 fhDZNormalized[i] = dynamic_cast<TH1F*> (gDirectory->Get("dZNormalized" ));
1379 fhDXYvsDZNormalized[i] = dynamic_cast<TH2F*> (gDirectory->Get("dXYvsDZNormalized"));
1380 fhNSigmaToVertex[i] = dynamic_cast<TH1F*> (gDirectory->Get("nSigmaToVertex"));
1382 fhPt[i] = dynamic_cast<TH1F*> (gDirectory->Get("pt"));
1383 fhEta[i] = dynamic_cast<TH1F*> (gDirectory->Get("eta"));
1385 gDirectory->cd("../");
1388 gDirectory->cd("..");
1393 //____________________________________________________________________
1394 void AliESDtrackCuts::SaveHistograms(const Char_t* dir) {
1396 // saves the histograms in a directory (dir)
1399 if (!fHistogramsOn) {
1400 AliDebug(0, "Histograms not on - cannot save histograms!!!");
1407 gDirectory->mkdir(dir);
1408 gDirectory->cd(dir);
1410 gDirectory->mkdir("before_cuts");
1411 gDirectory->mkdir("after_cuts");
1413 // a factor of 2 is needed since n sigma is positive
1414 ffDTheoretical->SetParameter(0,2*fhNSigmaToVertex[0]->Integral("width"));
1415 ffDTheoretical->Write("nSigmaToVertexTheory");
1417 fhCutStatistics->Write();
1418 fhCutCorrelation->Write();
1420 for (Int_t i=0; i<2; i++) {
1422 gDirectory->cd("before_cuts");
1424 gDirectory->cd("after_cuts");
1426 fhNClustersITS[i] ->Write();
1427 fhNClustersTPC[i] ->Write();
1428 fhChi2PerClusterITS[i] ->Write();
1429 fhChi2PerClusterTPC[i] ->Write();
1437 fhRel1PtUncertainty[i] ->Write();
1441 fhDXYDZ[i] ->Write();
1442 fhDXYvsDZ[i] ->Write();
1444 fhDXYNormalized[i] ->Write();
1445 fhDZNormalized[i] ->Write();
1446 fhDXYvsDZNormalized[i] ->Write();
1447 fhNSigmaToVertex[i] ->Write();
1452 gDirectory->cd("../");
1455 gDirectory->cd("../");
1458 //____________________________________________________________________
1459 void AliESDtrackCuts::DrawHistograms()
1461 // draws some histograms
1463 TCanvas* canvas1 = new TCanvas(Form("%s_1", GetName()), "Track Quality Results1", 800, 800);
1464 canvas1->Divide(2, 2);
1467 fhNClustersTPC[0]->SetStats(kFALSE);
1468 fhNClustersTPC[0]->Draw();
1471 fhChi2PerClusterTPC[0]->SetStats(kFALSE);
1472 fhChi2PerClusterTPC[0]->Draw();
1475 fhNSigmaToVertex[0]->SetStats(kFALSE);
1476 fhNSigmaToVertex[0]->GetXaxis()->SetRangeUser(0, 10);
1477 fhNSigmaToVertex[0]->Draw();
1479 canvas1->SaveAs(Form("%s_%s.gif", GetName(), canvas1->GetName()));
1481 TCanvas* canvas2 = new TCanvas(Form("%s_2", GetName()), "Track Quality Results2", 1200, 800);
1482 canvas2->Divide(3, 2);
1485 fhC11[0]->SetStats(kFALSE);
1490 fhC22[0]->SetStats(kFALSE);
1495 fhC33[0]->SetStats(kFALSE);
1500 fhC44[0]->SetStats(kFALSE);
1505 fhC55[0]->SetStats(kFALSE);
1510 fhRel1PtUncertainty[0]->SetStats(kFALSE);
1512 fhRel1PtUncertainty[0]->Draw();
1514 canvas2->SaveAs(Form("%s_%s.gif", GetName(), canvas2->GetName()));
1516 TCanvas* canvas3 = new TCanvas(Form("%s_3", GetName()), "Track Quality Results3", 1200, 800);
1517 canvas3->Divide(3, 2);
1520 fhDXY[0]->SetStats(kFALSE);
1525 fhDZ[0]->SetStats(kFALSE);
1530 fhDXYvsDZ[0]->SetStats(kFALSE);
1532 gPad->SetRightMargin(0.15);
1533 fhDXYvsDZ[0]->Draw("COLZ");
1536 fhDXYNormalized[0]->SetStats(kFALSE);
1538 fhDXYNormalized[0]->Draw();
1541 fhDZNormalized[0]->SetStats(kFALSE);
1543 fhDZNormalized[0]->Draw();
1546 fhDXYvsDZNormalized[0]->SetStats(kFALSE);
1548 gPad->SetRightMargin(0.15);
1549 fhDXYvsDZNormalized[0]->Draw("COLZ");
1551 canvas3->SaveAs(Form("%s_%s.gif", GetName(), canvas3->GetName()));
1553 TCanvas* canvas4 = new TCanvas(Form("%s_4", GetName()), "Track Quality Results4", 800, 500);
1554 canvas4->Divide(2, 1);
1557 fhCutStatistics->SetStats(kFALSE);
1558 fhCutStatistics->LabelsOption("v");
1559 gPad->SetBottomMargin(0.3);
1560 fhCutStatistics->Draw();
1563 fhCutCorrelation->SetStats(kFALSE);
1564 fhCutCorrelation->LabelsOption("v");
1565 gPad->SetBottomMargin(0.3);
1566 gPad->SetLeftMargin(0.3);
1567 fhCutCorrelation->Draw("COLZ");
1569 canvas4->SaveAs(Form("%s_%s.gif", GetName(), canvas4->GetName()));
1572 fhDXYvsDZNormalized[0]->SetStats(kFALSE);
1573 fhDXYvsDZNormalized[0]->DrawCopy("COLZ");
1576 fhNClustersTPC[0]->SetStats(kFALSE);
1577 fhNClustersTPC[0]->DrawCopy();
1580 fhChi2PerClusterITS[0]->SetStats(kFALSE);
1581 fhChi2PerClusterITS[0]->DrawCopy();
1582 fhChi2PerClusterITS[1]->SetLineColor(2);
1583 fhChi2PerClusterITS[1]->DrawCopy("SAME");
1586 fhChi2PerClusterTPC[0]->SetStats(kFALSE);
1587 fhChi2PerClusterTPC[0]->DrawCopy();
1588 fhChi2PerClusterTPC[1]->SetLineColor(2);
1589 fhChi2PerClusterTPC[1]->DrawCopy("SAME");*/
1591 //--------------------------------------------------------------------------
1592 void AliESDtrackCuts::SetPtDepDCACuts(Double_t pt) {
1594 // set the pt-dependent DCA cuts
1597 if(f1CutMaxDCAToVertexXYPtDep) {
1598 fCutMaxDCAToVertexXY=f1CutMaxDCAToVertexXYPtDep->Eval(pt);
1601 if(f1CutMaxDCAToVertexZPtDep) {
1602 fCutMaxDCAToVertexZ=f1CutMaxDCAToVertexZPtDep->Eval(pt);
1605 if(f1CutMinDCAToVertexXYPtDep) {
1606 fCutMinDCAToVertexXY=f1CutMinDCAToVertexXYPtDep->Eval(pt);
1609 if(f1CutMinDCAToVertexZPtDep) {
1610 fCutMinDCAToVertexZ=f1CutMinDCAToVertexZPtDep->Eval(pt);
1619 //--------------------------------------------------------------------------
1620 Bool_t AliESDtrackCuts::CheckPtDepDCA(TString dist,Bool_t print) const {
1622 // Check the correctness of the string syntax
1624 Bool_t retval=kTRUE;
1626 if(!dist.Contains("pt")) {
1627 if(print) printf("AliESDtrackCuts::CheckPtDepDCA(): string must contain \"pt\"\n");
1633 void AliESDtrackCuts::SetMaxDCAToVertexXYPtDep(const char *dist){
1635 if(f1CutMaxDCAToVertexXYPtDep){
1636 delete f1CutMaxDCAToVertexXYPtDep;
1638 f1CutMaxDCAToVertexXYPtDep = 0;
1639 fCutMaxDCAToVertexXYPtDep = "";
1641 if(!CheckPtDepDCA(dist,kTRUE)){
1644 fCutMaxDCAToVertexXYPtDep = dist;
1646 tmp.ReplaceAll("pt","x");
1647 f1CutMaxDCAToVertexXYPtDep = new TFormula("f1CutMaxDCAToVertexXYPtDep",tmp.Data());
1651 void AliESDtrackCuts::SetMaxDCAToVertexZPtDep(const char *dist){
1654 if(f1CutMaxDCAToVertexZPtDep){
1655 delete f1CutMaxDCAToVertexZPtDep;
1657 f1CutMaxDCAToVertexZPtDep = 0;
1658 fCutMaxDCAToVertexZPtDep = "";
1660 if(!CheckPtDepDCA(dist,kTRUE))return;
1662 fCutMaxDCAToVertexZPtDep = dist;
1664 tmp.ReplaceAll("pt","x");
1665 f1CutMaxDCAToVertexZPtDep = new TFormula("f1CutMaxDCAToVertexZPtDep",tmp.Data());
1671 void AliESDtrackCuts::SetMinDCAToVertexXYPtDep(const char *dist){
1674 if(f1CutMinDCAToVertexXYPtDep){
1675 delete f1CutMinDCAToVertexXYPtDep;
1677 f1CutMinDCAToVertexXYPtDep = 0;
1678 fCutMinDCAToVertexXYPtDep = "";
1680 if(!CheckPtDepDCA(dist,kTRUE))return;
1682 fCutMinDCAToVertexXYPtDep = dist;
1684 tmp.ReplaceAll("pt","x");
1685 f1CutMinDCAToVertexXYPtDep = new TFormula("f1CutMinDCAToVertexXYPtDep",tmp.Data());
1690 void AliESDtrackCuts::SetMinDCAToVertexZPtDep(const char *dist){
1694 if(f1CutMinDCAToVertexZPtDep){
1695 delete f1CutMinDCAToVertexZPtDep;
1697 f1CutMinDCAToVertexZPtDep = 0;
1698 fCutMinDCAToVertexZPtDep = "";
1700 if(!CheckPtDepDCA(dist,kTRUE))return;
1701 fCutMinDCAToVertexZPtDep = dist;
1703 tmp.ReplaceAll("pt","x");
1704 f1CutMinDCAToVertexZPtDep = new TFormula("f1CutMinDCAToVertexZPtDep",tmp.Data());