Line data Source code
1 : /**************************************************************************
2 : * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 : * *
4 : * Author: The ALICE Off-line Project. *
5 : * Contributors are mentioned in the code where appropriate. *
6 : * *
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 : **************************************************************************/
15 :
16 : /* $Id: AliESDtrackCuts.cxx 24534 2008-03-16 22:22:11Z fca $ */
17 :
18 : #include "AliESDtrackCuts.h"
19 :
20 : #include <AliESDtrack.h>
21 : #include <AliESDVertex.h>
22 : #include <AliESDEvent.h>
23 : #include <AliMultiplicity.h>
24 : #include <AliLog.h>
25 :
26 : #include <TTree.h>
27 : #include <TCanvas.h>
28 : #include <TDirectory.h>
29 : #include <TH2F.h>
30 : #include <TF1.h>
31 : #include <TBits.h>
32 :
33 : //____________________________________________________________________
34 170 : ClassImp(AliESDtrackCuts)
35 :
36 : // Cut names
37 : const Char_t* AliESDtrackCuts::fgkCutNames[kNCuts] = {
38 : "require TPC refit",
39 : "require TPC standalone",
40 : "require ITS refit",
41 : "n clusters TPC",
42 : "n clusters ITS",
43 : "#Chi^{2}/cluster TPC",
44 : "#Chi^{2}/cluster ITS",
45 : "cov 11",
46 : "cov 22",
47 : "cov 33",
48 : "cov 44",
49 : "cov 55",
50 : "trk-to-vtx",
51 : "trk-to-vtx failed",
52 : "kink daughters",
53 : "p",
54 : "p_{T}",
55 : "p_{x}",
56 : "p_{y}",
57 : "p_{z}",
58 : "eta",
59 : "y",
60 : "trk-to-vtx max dca 2D absolute",
61 : "trk-to-vtx max dca xy absolute",
62 : "trk-to-vtx max dca z absolute",
63 : "trk-to-vtx min dca 2D absolute",
64 : "trk-to-vtx min dca xy absolute",
65 : "trk-to-vtx min dca z absolute",
66 : "SPD cluster requirement",
67 : "SDD cluster requirement",
68 : "SSD cluster requirement",
69 : "require ITS stand-alone",
70 : "rel 1/pt uncertainty",
71 : "TPC n shared clusters",
72 : "TPC rel shared clusters",
73 : "require ITS Pid",
74 : "n crossed rows TPC",
75 : "n crossed rows / n findable clusters",
76 : "missing ITS points",
77 : "#Chi^{2} TPC constrained vs. global",
78 : "require TOF out",
79 : "TOF Distance cut",
80 : "min length in active volume TPC",
81 : "n-geometrical+n-crossed-row and n-clusters cut",
82 : "track is in distorted TPC region"
83 : };
84 :
85 : AliESDtrackCuts* AliESDtrackCuts::fgMultEstTrackCuts[AliESDtrackCuts::kNMultEstTrackCuts] = { 0, 0, 0, 0 };
86 : Char_t AliESDtrackCuts::fgBeamTypeFlag = -1;
87 :
88 : //____________________________________________________________________
89 : AliESDtrackCuts::AliESDtrackCuts(const Char_t* name, const Char_t* title) :
90 22 : AliAnalysisCuts(name,title),
91 22 : fCutMinNClusterTPC(0),
92 22 : fCutMinNClusterITS(0),
93 22 : fCutMinNCrossedRowsTPC(0),
94 22 : fCutMinRatioCrossedRowsOverFindableClustersTPC(0),
95 22 : f1CutMinNClustersTPCPtDep(0x0),
96 22 : fCutMaxPtDepNClustersTPC(0),
97 22 : fCutMinLengthActiveVolumeTPC(0),
98 22 : fDeadZoneWidth(0), // width of the TPC dead zone (missing pads + PRF +ExB)
99 22 : fCutGeoNcrNclLength(0), // cut on the geometical length condition Ngeom(cm)>cutGeoNcrNclLength default=130
100 22 : fCutGeoNcrNclGeom1Pt(0), // 1/pt dependence slope cutGeoNcrNclLength:=fCutGeoNcrNclLength-abs(1/pt)^fCutGeoNcrNclGeom1Pt
101 22 : fCutGeoNcrNclFractionNcr(0), // relative fraction cut Ncr condition Ncr>cutGeoNcrNclFractionNcr*fCutGeoNcrNclLength
102 22 : fCutGeoNcrNclFractionNcl(0), // relative fraction cut Ncr condition Ncl>cutGeoNcrNclFractionNcl
103 22 : fCutOutDistortedRegionTPC(kFALSE), // flag if distorted TPC regions should be cut out
104 22 : fCutMaxChi2PerClusterTPC(0),
105 22 : fCutMaxChi2PerClusterITS(0),
106 22 : fCutMaxChi2TPCConstrainedVsGlobal(0),
107 22 : fCutMaxChi2TPCConstrainedVsGlobalVertexType(kVertexTracks | kVertexSPD),
108 22 : fCutMaxMissingITSPoints(0),
109 22 : fCutMaxC11(0),
110 22 : fCutMaxC22(0),
111 22 : fCutMaxC33(0),
112 22 : fCutMaxC44(0),
113 22 : fCutMaxC55(0),
114 22 : fCutMaxRel1PtUncertainty(0),
115 22 : fCutAcceptKinkDaughters(0),
116 22 : fCutAcceptSharedTPCClusters(0),
117 22 : fCutMaxFractionSharedTPCClusters(0),
118 22 : fCutRequireTPCRefit(0),
119 22 : fCutRequireTPCStandAlone(0),
120 22 : fCutRequireITSRefit(0),
121 22 : fCutRequireITSPid(0),
122 22 : fCutRequireITSStandAlone(0),
123 22 : fCutRequireITSpureSA(0),
124 22 : fCutNsigmaToVertex(0),
125 22 : fCutSigmaToVertexRequired(0),
126 22 : fCutMaxDCAToVertexXY(0),
127 22 : fCutMaxDCAToVertexZ(0),
128 22 : fCutMinDCAToVertexXY(0),
129 22 : fCutMinDCAToVertexZ(0),
130 22 : fCutMaxDCAToVertexXYPtDep(""),
131 22 : fCutMaxDCAToVertexZPtDep(""),
132 22 : fCutMinDCAToVertexXYPtDep(""),
133 22 : fCutMinDCAToVertexZPtDep(""),
134 22 : f1CutMaxDCAToVertexXYPtDep(0x0),
135 22 : f1CutMaxDCAToVertexZPtDep(0x0),
136 22 : f1CutMinDCAToVertexXYPtDep(0x0),
137 22 : f1CutMinDCAToVertexZPtDep(0x0),
138 22 : fCutDCAToVertex2D(0),
139 22 : fPMin(0),
140 22 : fPMax(0),
141 22 : fPtMin(0),
142 22 : fPtMax(0),
143 22 : fPxMin(0),
144 22 : fPxMax(0),
145 22 : fPyMin(0),
146 22 : fPyMax(0),
147 22 : fPzMin(0),
148 22 : fPzMax(0),
149 22 : fEtaMin(0),
150 22 : fEtaMax(0),
151 22 : fRapMin(0),
152 22 : fRapMax(0),
153 22 : fCutRequireTOFout(kFALSE),
154 22 : fFlagCutTOFdistance(kFALSE),
155 22 : fCutTOFdistance(3.),
156 22 : fHistogramsOn(0),
157 22 : ffDTheoretical(0),
158 22 : fhCutStatistics(0),
159 22 : fhCutCorrelation(0)
160 98 : {
161 : //
162 : // constructor
163 : //
164 :
165 22 : Init();
166 :
167 : //##############################################################################
168 : // setting default cuts
169 22 : SetMinNClustersTPC();
170 22 : SetMinNClustersITS();
171 22 : SetMinNCrossedRowsTPC();
172 22 : SetMinRatioCrossedRowsOverFindableClustersTPC();
173 22 : SetMaxChi2PerClusterTPC();
174 22 : SetMaxChi2PerClusterITS();
175 22 : SetMaxChi2TPCConstrainedGlobal();
176 22 : SetMaxChi2TPCConstrainedGlobalVertexType();
177 22 : SetMaxNOfMissingITSPoints();
178 22 : SetMaxCovDiagonalElements();
179 22 : SetMaxRel1PtUncertainty();
180 22 : SetRequireTPCRefit();
181 22 : SetRequireTPCStandAlone();
182 22 : SetRequireITSRefit();
183 22 : SetRequireITSPid(kFALSE);
184 22 : SetRequireITSStandAlone(kFALSE);
185 22 : SetRequireITSPureStandAlone(kFALSE);
186 22 : SetAcceptKinkDaughters();
187 22 : SetAcceptSharedTPCClusters();
188 22 : SetMaxFractionSharedTPCClusters();
189 22 : SetMaxNsigmaToVertex();
190 22 : SetMaxDCAToVertexXY();
191 22 : SetMaxDCAToVertexZ();
192 22 : SetDCAToVertex2D();
193 22 : SetMinDCAToVertexXY();
194 22 : SetMinDCAToVertexZ();
195 22 : SetPRange();
196 22 : SetPtRange();
197 22 : SetPxRange();
198 22 : SetPyRange();
199 22 : SetPzRange();
200 22 : SetEtaRange();
201 22 : SetRapRange();
202 22 : SetClusterRequirementITS(kSPD);
203 22 : SetClusterRequirementITS(kSDD);
204 22 : SetClusterRequirementITS(kSSD);
205 :
206 22 : SetHistogramsOn();
207 38 : }
208 :
209 : //_____________________________________________________________________________
210 : AliESDtrackCuts::AliESDtrackCuts(const AliESDtrackCuts &c) :
211 0 : AliAnalysisCuts(c),
212 0 : fCutMinNClusterTPC(0),
213 0 : fCutMinNClusterITS(0),
214 0 : fCutMinNCrossedRowsTPC(0),
215 0 : fCutMinRatioCrossedRowsOverFindableClustersTPC(0),
216 0 : f1CutMinNClustersTPCPtDep(0x0),
217 0 : fCutMaxPtDepNClustersTPC(0),
218 0 : fCutMinLengthActiveVolumeTPC(0),
219 0 : fDeadZoneWidth(0), // width of the TPC dead zone (missing pads + PRF +ExB)
220 0 : fCutGeoNcrNclLength(0), // cut on the geometical length condition Ngeom(cm)>cutGeoNcrNclLength default=130
221 0 : fCutGeoNcrNclGeom1Pt(0), // 1/pt dependence slope cutGeoNcrNclLength:=fCutGeoNcrNclLength-abs(1/pt)^fCutGeoNcrNclGeom1Pt
222 0 : fCutGeoNcrNclFractionNcr(0), // relative fraction cut Ncr condition Ncr>cutGeoNcrNclFractionNcr*fCutGeoNcrNclLength
223 0 : fCutGeoNcrNclFractionNcl(0), // relative fraction cut Ncr condition Ncl>cutGeoNcrNclFractionNcl
224 0 : fCutOutDistortedRegionTPC(kFALSE),
225 : //
226 0 : fCutMaxChi2PerClusterTPC(0),
227 0 : fCutMaxChi2PerClusterITS(0),
228 0 : fCutMaxChi2TPCConstrainedVsGlobal(0),
229 0 : fCutMaxChi2TPCConstrainedVsGlobalVertexType(kVertexTracks | kVertexSPD),
230 0 : fCutMaxMissingITSPoints(0),
231 0 : fCutMaxC11(0),
232 0 : fCutMaxC22(0),
233 0 : fCutMaxC33(0),
234 0 : fCutMaxC44(0),
235 0 : fCutMaxC55(0),
236 0 : fCutMaxRel1PtUncertainty(0),
237 0 : fCutAcceptKinkDaughters(0),
238 0 : fCutAcceptSharedTPCClusters(0),
239 0 : fCutMaxFractionSharedTPCClusters(0),
240 0 : fCutRequireTPCRefit(0),
241 0 : fCutRequireTPCStandAlone(0),
242 0 : fCutRequireITSRefit(0),
243 0 : fCutRequireITSPid(0),
244 0 : fCutRequireITSStandAlone(0),
245 0 : fCutRequireITSpureSA(0),
246 0 : fCutNsigmaToVertex(0),
247 0 : fCutSigmaToVertexRequired(0),
248 0 : fCutMaxDCAToVertexXY(0),
249 0 : fCutMaxDCAToVertexZ(0),
250 0 : fCutMinDCAToVertexXY(0),
251 0 : fCutMinDCAToVertexZ(0),
252 0 : fCutMaxDCAToVertexXYPtDep(""),
253 0 : fCutMaxDCAToVertexZPtDep(""),
254 0 : fCutMinDCAToVertexXYPtDep(""),
255 0 : fCutMinDCAToVertexZPtDep(""),
256 0 : f1CutMaxDCAToVertexXYPtDep(0x0),
257 0 : f1CutMaxDCAToVertexZPtDep(0x0),
258 0 : f1CutMinDCAToVertexXYPtDep(0x0),
259 0 : f1CutMinDCAToVertexZPtDep(0x0),
260 0 : fCutDCAToVertex2D(0),
261 0 : fPMin(0),
262 0 : fPMax(0),
263 0 : fPtMin(0),
264 0 : fPtMax(0),
265 0 : fPxMin(0),
266 0 : fPxMax(0),
267 0 : fPyMin(0),
268 0 : fPyMax(0),
269 0 : fPzMin(0),
270 0 : fPzMax(0),
271 0 : fEtaMin(0),
272 0 : fEtaMax(0),
273 0 : fRapMin(0),
274 0 : fRapMax(0),
275 0 : fCutRequireTOFout(kFALSE),
276 0 : fFlagCutTOFdistance(kFALSE),
277 0 : fCutTOFdistance(3.),
278 0 : fHistogramsOn(0),
279 0 : ffDTheoretical(0),
280 0 : fhCutStatistics(0),
281 0 : fhCutCorrelation(0)
282 0 : {
283 : //
284 : // copy constructor
285 : //
286 :
287 0 : ((AliESDtrackCuts &) c).Copy(*this);
288 0 : }
289 :
290 : AliESDtrackCuts::~AliESDtrackCuts()
291 24 : {
292 : //
293 : // destructor
294 : //
295 :
296 24 : for (Int_t i=0; i<2; i++) {
297 :
298 8 : if (fhNClustersITS[i])
299 16 : delete fhNClustersITS[i];
300 8 : if (fhNClustersTPC[i])
301 16 : delete fhNClustersTPC[i];
302 8 : if (fhNSharedClustersTPC[i])
303 16 : delete fhNSharedClustersTPC[i];
304 8 : if (fhNCrossedRowsTPC[i])
305 16 : delete fhNCrossedRowsTPC[i];
306 8 : if (fhRatioCrossedRowsOverFindableClustersTPC[i])
307 16 : delete fhRatioCrossedRowsOverFindableClustersTPC[i];
308 8 : if (fhChi2PerClusterITS[i])
309 16 : delete fhChi2PerClusterITS[i];
310 8 : if (fhChi2PerClusterTPC[i])
311 16 : delete fhChi2PerClusterTPC[i];
312 8 : if (fhChi2TPCConstrainedVsGlobal[i])
313 16 : delete fhChi2TPCConstrainedVsGlobal[i];
314 8 : if(fhNClustersForITSPID[i])
315 16 : delete fhNClustersForITSPID[i];
316 8 : if(fhNMissingITSPoints[i])
317 16 : delete fhNMissingITSPoints[i];
318 8 : if (fhC11[i])
319 16 : delete fhC11[i];
320 8 : if (fhC22[i])
321 16 : delete fhC22[i];
322 8 : if (fhC33[i])
323 16 : delete fhC33[i];
324 8 : if (fhC44[i])
325 16 : delete fhC44[i];
326 8 : if (fhC55[i])
327 16 : delete fhC55[i];
328 :
329 8 : if (fhRel1PtUncertainty[i])
330 16 : delete fhRel1PtUncertainty[i];
331 :
332 8 : if (fhDXY[i])
333 16 : delete fhDXY[i];
334 8 : if (fhDZ[i])
335 16 : delete fhDZ[i];
336 8 : if (fhDXYDZ[i])
337 16 : delete fhDXYDZ[i];
338 8 : if (fhDXYvsDZ[i])
339 16 : delete fhDXYvsDZ[i];
340 :
341 8 : if (fhDXYNormalized[i])
342 16 : delete fhDXYNormalized[i];
343 8 : if (fhDZNormalized[i])
344 16 : delete fhDZNormalized[i];
345 8 : if (fhDXYvsDZNormalized[i])
346 16 : delete fhDXYvsDZNormalized[i];
347 8 : if (fhNSigmaToVertex[i])
348 16 : delete fhNSigmaToVertex[i];
349 8 : if (fhPt[i])
350 16 : delete fhPt[i];
351 8 : if (fhEta[i])
352 16 : delete fhEta[i];
353 8 : if (fhTOFdistance[i])
354 16 : delete fhTOFdistance[i];
355 : }
356 :
357 4 : if(f1CutMaxDCAToVertexXYPtDep)delete f1CutMaxDCAToVertexXYPtDep;
358 4 : f1CutMaxDCAToVertexXYPtDep = 0;
359 4 : if( f1CutMaxDCAToVertexZPtDep) delete f1CutMaxDCAToVertexZPtDep;
360 4 : f1CutMaxDCAToVertexZPtDep = 0;
361 4 : if( f1CutMinDCAToVertexXYPtDep)delete f1CutMinDCAToVertexXYPtDep;
362 4 : f1CutMinDCAToVertexXYPtDep = 0;
363 4 : if(f1CutMinDCAToVertexZPtDep)delete f1CutMinDCAToVertexZPtDep;
364 4 : f1CutMinDCAToVertexZPtDep = 0;
365 :
366 :
367 4 : if (ffDTheoretical)
368 8 : delete ffDTheoretical;
369 :
370 4 : if (fhCutStatistics)
371 8 : delete fhCutStatistics;
372 4 : if (fhCutCorrelation)
373 8 : delete fhCutCorrelation;
374 :
375 4 : if(f1CutMinNClustersTPCPtDep)
376 0 : delete f1CutMinNClustersTPCPtDep;
377 :
378 12 : }
379 :
380 : void AliESDtrackCuts::Init()
381 : {
382 : //
383 : // sets everything to zero
384 : //
385 :
386 44 : fCutMinNClusterTPC = 0;
387 22 : fCutMinNClusterITS = 0;
388 :
389 22 : fCutMaxChi2PerClusterTPC = 0;
390 22 : fCutMaxChi2PerClusterITS = 0;
391 22 : fCutMaxChi2TPCConstrainedVsGlobal = 0;
392 22 : fCutMaxChi2TPCConstrainedVsGlobalVertexType = kVertexTracks | kVertexSPD;
393 22 : fCutMaxMissingITSPoints = 0;
394 :
395 176 : for (Int_t i = 0; i < 3; i++)
396 66 : fCutClusterRequirementITS[i] = kOff;
397 :
398 22 : fCutMaxC11 = 0;
399 22 : fCutMaxC22 = 0;
400 22 : fCutMaxC33 = 0;
401 22 : fCutMaxC44 = 0;
402 22 : fCutMaxC55 = 0;
403 :
404 22 : fCutMaxRel1PtUncertainty = 0;
405 :
406 22 : fCutAcceptKinkDaughters = 0;
407 22 : fCutAcceptSharedTPCClusters = 0;
408 22 : fCutMaxFractionSharedTPCClusters = 0;
409 22 : fCutRequireTPCRefit = 0;
410 22 : fCutRequireTPCStandAlone = 0;
411 22 : fCutRequireITSRefit = 0;
412 22 : fCutRequireITSPid = 0;
413 22 : fCutRequireITSStandAlone = 0;
414 22 : fCutRequireITSpureSA = 0;
415 :
416 22 : fCutNsigmaToVertex = 0;
417 22 : fCutSigmaToVertexRequired = 0;
418 22 : fCutMaxDCAToVertexXY = 0;
419 22 : fCutMaxDCAToVertexZ = 0;
420 22 : fCutDCAToVertex2D = 0;
421 22 : fCutMinDCAToVertexXY = 0;
422 22 : fCutMinDCAToVertexZ = 0;
423 22 : fCutMaxDCAToVertexXYPtDep = "";
424 22 : fCutMaxDCAToVertexZPtDep = "";
425 22 : fCutMinDCAToVertexXYPtDep = "";
426 22 : fCutMinDCAToVertexZPtDep = "";
427 :
428 22 : if(f1CutMaxDCAToVertexXYPtDep)delete f1CutMaxDCAToVertexXYPtDep;
429 22 : f1CutMaxDCAToVertexXYPtDep = 0;
430 22 : if( f1CutMaxDCAToVertexXYPtDep) delete f1CutMaxDCAToVertexXYPtDep;
431 22 : f1CutMaxDCAToVertexXYPtDep = 0;
432 22 : if( f1CutMaxDCAToVertexZPtDep) delete f1CutMaxDCAToVertexZPtDep;
433 22 : f1CutMaxDCAToVertexZPtDep = 0;
434 22 : if( f1CutMinDCAToVertexXYPtDep)delete f1CutMinDCAToVertexXYPtDep;
435 22 : f1CutMinDCAToVertexXYPtDep = 0;
436 22 : if(f1CutMinDCAToVertexZPtDep)delete f1CutMinDCAToVertexZPtDep;
437 22 : f1CutMinDCAToVertexZPtDep = 0;
438 :
439 :
440 22 : fPMin = 0;
441 22 : fPMax = 0;
442 22 : fPtMin = 0;
443 22 : fPtMax = 0;
444 22 : fPxMin = 0;
445 22 : fPxMax = 0;
446 22 : fPyMin = 0;
447 22 : fPyMax = 0;
448 22 : fPzMin = 0;
449 22 : fPzMax = 0;
450 22 : fEtaMin = 0;
451 22 : fEtaMax = 0;
452 22 : fRapMin = 0;
453 22 : fRapMax = 0;
454 :
455 22 : fHistogramsOn = kFALSE;
456 :
457 132 : for (Int_t i=0; i<2; ++i)
458 : {
459 44 : fhNClustersITS[i] = 0;
460 44 : fhNClustersTPC[i] = 0;
461 44 : fhNSharedClustersTPC[i] = 0;
462 44 : fhNCrossedRowsTPC[i] = 0;
463 44 : fhRatioCrossedRowsOverFindableClustersTPC[i] = 0;
464 :
465 44 : fhChi2PerClusterITS[i] = 0;
466 44 : fhChi2PerClusterTPC[i] = 0;
467 44 : fhChi2TPCConstrainedVsGlobal[i] = 0;
468 44 : fhNClustersForITSPID[i] = 0;
469 44 : fhNMissingITSPoints[i] = 0;
470 :
471 44 : fhC11[i] = 0;
472 44 : fhC22[i] = 0;
473 44 : fhC33[i] = 0;
474 44 : fhC44[i] = 0;
475 44 : fhC55[i] = 0;
476 :
477 44 : fhRel1PtUncertainty[i] = 0;
478 :
479 44 : fhDXY[i] = 0;
480 44 : fhDZ[i] = 0;
481 44 : fhDXYDZ[i] = 0;
482 44 : fhDXYvsDZ[i] = 0;
483 :
484 44 : fhDXYNormalized[i] = 0;
485 44 : fhDZNormalized[i] = 0;
486 44 : fhDXYvsDZNormalized[i] = 0;
487 44 : fhNSigmaToVertex[i] = 0;
488 :
489 44 : fhPt[i] = 0;
490 44 : fhEta[i] = 0;
491 44 : fhTOFdistance[i] = 0;
492 : }
493 22 : ffDTheoretical = 0;
494 :
495 22 : fhCutStatistics = 0;
496 22 : fhCutCorrelation = 0;
497 22 : }
498 :
499 : //_____________________________________________________________________________
500 : AliESDtrackCuts &AliESDtrackCuts::operator=(const AliESDtrackCuts &c)
501 : {
502 : //
503 : // Assignment operator
504 : //
505 :
506 0 : if (this != &c) ((AliESDtrackCuts &) c).Copy(*this);
507 0 : return *this;
508 : }
509 :
510 : //_____________________________________________________________________________
511 : void AliESDtrackCuts::Copy(TObject &c) const
512 : {
513 : //
514 : // Copy function
515 : //
516 :
517 0 : AliESDtrackCuts& target = (AliESDtrackCuts &) c;
518 :
519 0 : target.Init();
520 :
521 0 : target.fCutMinNClusterTPC = fCutMinNClusterTPC;
522 0 : target.fCutMinNClusterITS = fCutMinNClusterITS;
523 0 : target.fCutMinNCrossedRowsTPC = fCutMinNCrossedRowsTPC;
524 0 : target.fCutMinRatioCrossedRowsOverFindableClustersTPC = fCutMinRatioCrossedRowsOverFindableClustersTPC;
525 0 : if(f1CutMinNClustersTPCPtDep){
526 0 : target.f1CutMinNClustersTPCPtDep = (TFormula*) f1CutMinNClustersTPCPtDep->Clone("f1CutMinNClustersTPCPtDep");
527 0 : }
528 0 : target.fCutMaxPtDepNClustersTPC = fCutMaxPtDepNClustersTPC;
529 0 : target.fCutMinLengthActiveVolumeTPC = fCutMinLengthActiveVolumeTPC;
530 0 : target.fDeadZoneWidth=fDeadZoneWidth; // width of the TPC dead zone (missing pads + PRF +ExB)
531 0 : target.fCutGeoNcrNclLength=fCutGeoNcrNclLength; // cut on the geometical length condition Ngeom(cm)>cutGeoNcrNclLength default=130
532 0 : target.fCutGeoNcrNclGeom1Pt=fCutGeoNcrNclGeom1Pt; // 1/pt dependence slope cutGeoNcrNclLength:=fCutGeoNcrNclLength-abs(1/pt)^fCutGeoNcrNclGeom1Pt
533 0 : target.fCutGeoNcrNclFractionNcr=fCutGeoNcrNclFractionNcr; // relative fraction cut Ncr condition Ncr>cutGeoNcrNclFractionNcr*fCutGeoNcrNclLength
534 0 : target.fCutGeoNcrNclFractionNcl=fCutGeoNcrNclFractionNcl; // relative fraction cut Ncr condition Ncl>cutGeoNcrNclFractionNcl
535 0 : target.fCutOutDistortedRegionTPC=fCutOutDistortedRegionTPC;
536 :
537 :
538 0 : target.fCutMaxChi2PerClusterTPC = fCutMaxChi2PerClusterTPC;
539 0 : target.fCutMaxChi2PerClusterITS = fCutMaxChi2PerClusterITS;
540 0 : target.fCutMaxChi2TPCConstrainedVsGlobal = fCutMaxChi2TPCConstrainedVsGlobal;
541 0 : target.fCutMaxChi2TPCConstrainedVsGlobalVertexType = fCutMaxChi2TPCConstrainedVsGlobalVertexType;
542 0 : target.fCutMaxMissingITSPoints = fCutMaxMissingITSPoints;
543 :
544 0 : for (Int_t i = 0; i < 3; i++)
545 0 : target.fCutClusterRequirementITS[i] = fCutClusterRequirementITS[i];
546 :
547 0 : target.fCutMaxC11 = fCutMaxC11;
548 0 : target.fCutMaxC22 = fCutMaxC22;
549 0 : target.fCutMaxC33 = fCutMaxC33;
550 0 : target.fCutMaxC44 = fCutMaxC44;
551 0 : target.fCutMaxC55 = fCutMaxC55;
552 :
553 0 : target.fCutMaxRel1PtUncertainty = fCutMaxRel1PtUncertainty;
554 :
555 0 : target.fCutAcceptKinkDaughters = fCutAcceptKinkDaughters;
556 0 : target.fCutAcceptSharedTPCClusters = fCutAcceptSharedTPCClusters;
557 0 : target.fCutMaxFractionSharedTPCClusters = fCutMaxFractionSharedTPCClusters;
558 0 : target.fCutRequireTPCRefit = fCutRequireTPCRefit;
559 0 : target.fCutRequireTPCStandAlone = fCutRequireTPCStandAlone;
560 0 : target.fCutRequireITSRefit = fCutRequireITSRefit;
561 0 : target.fCutRequireITSPid = fCutRequireITSPid;
562 0 : target.fCutRequireITSStandAlone = fCutRequireITSStandAlone;
563 0 : target.fCutRequireITSpureSA = fCutRequireITSpureSA;
564 :
565 0 : target.fCutNsigmaToVertex = fCutNsigmaToVertex;
566 0 : target.fCutSigmaToVertexRequired = fCutSigmaToVertexRequired;
567 0 : target.fCutMaxDCAToVertexXY = fCutMaxDCAToVertexXY;
568 0 : target.fCutMaxDCAToVertexZ = fCutMaxDCAToVertexZ;
569 0 : target.fCutDCAToVertex2D = fCutDCAToVertex2D;
570 0 : target.fCutMinDCAToVertexXY = fCutMinDCAToVertexXY;
571 0 : target.fCutMinDCAToVertexZ = fCutMinDCAToVertexZ;
572 :
573 0 : target.fCutMaxDCAToVertexXYPtDep = fCutMaxDCAToVertexXYPtDep;
574 0 : if(fCutMaxDCAToVertexXYPtDep.Length()>0)target.SetMaxDCAToVertexXYPtDep(fCutMaxDCAToVertexXYPtDep.Data());
575 :
576 0 : target.fCutMaxDCAToVertexZPtDep = fCutMaxDCAToVertexZPtDep;
577 0 : if(fCutMaxDCAToVertexZPtDep.Length()>0)target.SetMaxDCAToVertexZPtDep(fCutMaxDCAToVertexZPtDep.Data());
578 :
579 0 : target.fCutMinDCAToVertexXYPtDep = fCutMinDCAToVertexXYPtDep;
580 0 : if(fCutMinDCAToVertexXYPtDep.Length()>0)target.SetMinDCAToVertexXYPtDep(fCutMinDCAToVertexXYPtDep.Data());
581 :
582 0 : target.fCutMinDCAToVertexZPtDep = fCutMinDCAToVertexZPtDep;
583 0 : if(fCutMinDCAToVertexZPtDep.Length()>0)target.SetMinDCAToVertexZPtDep(fCutMinDCAToVertexZPtDep.Data());
584 :
585 0 : target.fPMin = fPMin;
586 0 : target.fPMax = fPMax;
587 0 : target.fPtMin = fPtMin;
588 0 : target.fPtMax = fPtMax;
589 0 : target.fPxMin = fPxMin;
590 0 : target.fPxMax = fPxMax;
591 0 : target.fPyMin = fPyMin;
592 0 : target.fPyMax = fPyMax;
593 0 : target.fPzMin = fPzMin;
594 0 : target.fPzMax = fPzMax;
595 0 : target.fEtaMin = fEtaMin;
596 0 : target.fEtaMax = fEtaMax;
597 0 : target.fRapMin = fRapMin;
598 0 : target.fRapMax = fRapMax;
599 :
600 0 : target.fFlagCutTOFdistance = fFlagCutTOFdistance;
601 0 : target.fCutTOFdistance = fCutTOFdistance;
602 0 : target.fCutRequireTOFout = fCutRequireTOFout;
603 :
604 0 : target.fHistogramsOn = fHistogramsOn;
605 :
606 0 : for (Int_t i=0; i<2; ++i)
607 : {
608 0 : if (fhNClustersITS[i]) target.fhNClustersITS[i] = (TH1F*) fhNClustersITS[i]->Clone();
609 0 : if (fhNClustersTPC[i]) target.fhNClustersTPC[i] = (TH1F*) fhNClustersTPC[i]->Clone();
610 0 : if (fhNSharedClustersTPC[i]) target.fhNSharedClustersTPC[i] = (TH1F*) fhNSharedClustersTPC[i]->Clone();
611 0 : if (fhNCrossedRowsTPC[i]) target.fhNCrossedRowsTPC[i] = (TH1F*) fhNCrossedRowsTPC[i]->Clone();
612 0 : if (fhRatioCrossedRowsOverFindableClustersTPC[i]) target.fhRatioCrossedRowsOverFindableClustersTPC[i] = (TH1F*) fhRatioCrossedRowsOverFindableClustersTPC[i]->Clone();
613 :
614 0 : if (fhChi2PerClusterITS[i]) target.fhChi2PerClusterITS[i] = (TH1F*) fhChi2PerClusterITS[i]->Clone();
615 0 : if (fhChi2PerClusterTPC[i]) target.fhChi2PerClusterTPC[i] = (TH1F*) fhChi2PerClusterTPC[i]->Clone();
616 0 : if (fhChi2TPCConstrainedVsGlobal[i]) target.fhChi2TPCConstrainedVsGlobal[i] = (TH1F*) fhChi2TPCConstrainedVsGlobal[i]->Clone();
617 0 : if (fhNClustersForITSPID[i]) target.fhNClustersForITSPID[i] = (TH1F*) fhNClustersForITSPID[i]->Clone();
618 0 : if (fhNMissingITSPoints[i]) target.fhNMissingITSPoints[i] = (TH1F*) fhNMissingITSPoints[i]->Clone();
619 :
620 0 : if (fhC11[i]) target.fhC11[i] = (TH1F*) fhC11[i]->Clone();
621 0 : if (fhC22[i]) target.fhC22[i] = (TH1F*) fhC22[i]->Clone();
622 0 : if (fhC33[i]) target.fhC33[i] = (TH1F*) fhC33[i]->Clone();
623 0 : if (fhC44[i]) target.fhC44[i] = (TH1F*) fhC44[i]->Clone();
624 0 : if (fhC55[i]) target.fhC55[i] = (TH1F*) fhC55[i]->Clone();
625 :
626 0 : if (fhRel1PtUncertainty[i]) target.fhRel1PtUncertainty[i] = (TH1F*) fhRel1PtUncertainty[i]->Clone();
627 :
628 0 : if (fhDXY[i]) target.fhDXY[i] = (TH1F*) fhDXY[i]->Clone();
629 0 : if (fhDZ[i]) target.fhDZ[i] = (TH1F*) fhDZ[i]->Clone();
630 0 : if (fhDXYDZ[i]) target.fhDXYDZ[i] = (TH1F*) fhDXYDZ[i]->Clone();
631 0 : if (fhDXYvsDZ[i]) target.fhDXYvsDZ[i] = (TH2F*) fhDXYvsDZ[i]->Clone();
632 :
633 0 : if (fhDXYNormalized[i]) target.fhDXYNormalized[i] = (TH1F*) fhDXYNormalized[i]->Clone();
634 0 : if (fhDZNormalized[i]) target.fhDZNormalized[i] = (TH1F*) fhDZNormalized[i]->Clone();
635 0 : if (fhDXYvsDZNormalized[i]) target.fhDXYvsDZNormalized[i] = (TH2F*) fhDXYvsDZNormalized[i]->Clone();
636 0 : if (fhNSigmaToVertex[i]) target.fhNSigmaToVertex[i] = (TH1F*) fhNSigmaToVertex[i]->Clone();
637 :
638 0 : if (fhPt[i]) target.fhPt[i] = (TH1F*) fhPt[i]->Clone();
639 0 : if (fhEta[i]) target.fhEta[i] = (TH1F*) fhEta[i]->Clone();
640 0 : if (fhTOFdistance[i]) target.fhTOFdistance[i] = (TH2F*) fhTOFdistance[i]->Clone();
641 : }
642 0 : if (ffDTheoretical) target.ffDTheoretical = (TF1*) ffDTheoretical->Clone();
643 :
644 0 : if (fhCutStatistics) target.fhCutStatistics = (TH1F*) fhCutStatistics->Clone();
645 0 : if (fhCutCorrelation) target.fhCutCorrelation = (TH2F*) fhCutCorrelation->Clone();
646 :
647 0 : TNamed::Copy(c);
648 0 : }
649 :
650 : //_____________________________________________________________________________
651 : Long64_t AliESDtrackCuts::Merge(TCollection* list) {
652 : // Merge a list of AliESDtrackCuts objects with this (needed for PROOF)
653 : // Returns the number of merged objects (including this)
654 0 : if (!list)
655 0 : return 0;
656 0 : if (list->IsEmpty())
657 0 : return 1;
658 0 : if (!fHistogramsOn)
659 0 : return 0;
660 0 : TIterator* iter = list->MakeIterator();
661 : TObject* obj;
662 :
663 : // collection of measured and generated histograms
664 : Int_t count = 0;
665 0 : while ((obj = iter->Next())) {
666 :
667 0 : AliESDtrackCuts* entry = dynamic_cast<AliESDtrackCuts*>(obj);
668 0 : if (entry == 0)
669 0 : continue;
670 :
671 0 : if (!entry->fHistogramsOn)
672 0 : continue;
673 :
674 0 : for (Int_t i=0; i<2; i++) {
675 :
676 0 : fhNClustersITS[i] ->Add(entry->fhNClustersITS[i] );
677 0 : fhNClustersTPC[i] ->Add(entry->fhNClustersTPC[i] );
678 0 : if (fhNSharedClustersTPC[i])
679 0 : fhNSharedClustersTPC[i] ->Add(entry->fhNSharedClustersTPC[i] );
680 0 : if (fhNCrossedRowsTPC[i])
681 0 : fhNCrossedRowsTPC[i] ->Add(entry->fhNCrossedRowsTPC[i] );
682 0 : if (fhRatioCrossedRowsOverFindableClustersTPC[i])
683 0 : fhRatioCrossedRowsOverFindableClustersTPC[i] ->Add(entry->fhRatioCrossedRowsOverFindableClustersTPC[i] );
684 :
685 0 : fhChi2PerClusterITS[i] ->Add(entry->fhChi2PerClusterITS[i]);
686 0 : fhChi2PerClusterTPC[i] ->Add(entry->fhChi2PerClusterTPC[i]);
687 0 : if (fhChi2TPCConstrainedVsGlobal[i])
688 0 : fhChi2TPCConstrainedVsGlobal[i]->Add(entry->fhChi2TPCConstrainedVsGlobal[i]);
689 0 : if (fhNClustersForITSPID[i])
690 0 : fhNClustersForITSPID[i]->Add(entry->fhNClustersForITSPID[i]);
691 0 : if (fhNMissingITSPoints[i])
692 0 : fhNMissingITSPoints[i] ->Add(entry->fhNMissingITSPoints[i]);
693 :
694 0 : fhC11[i] ->Add(entry->fhC11[i] );
695 0 : fhC22[i] ->Add(entry->fhC22[i] );
696 0 : fhC33[i] ->Add(entry->fhC33[i] );
697 0 : fhC44[i] ->Add(entry->fhC44[i] );
698 0 : fhC55[i] ->Add(entry->fhC55[i] );
699 :
700 0 : fhRel1PtUncertainty[i] ->Add(entry->fhRel1PtUncertainty[i]);
701 :
702 0 : fhDXY[i] ->Add(entry->fhDXY[i] );
703 0 : fhDZ[i] ->Add(entry->fhDZ[i] );
704 0 : fhDXYDZ[i] ->Add(entry->fhDXYDZ[i] );
705 0 : fhDXYvsDZ[i] ->Add(entry->fhDXYvsDZ[i] );
706 :
707 0 : fhDXYNormalized[i] ->Add(entry->fhDXYNormalized[i] );
708 0 : fhDZNormalized[i] ->Add(entry->fhDZNormalized[i] );
709 0 : fhDXYvsDZNormalized[i] ->Add(entry->fhDXYvsDZNormalized[i]);
710 0 : fhNSigmaToVertex[i] ->Add(entry->fhNSigmaToVertex[i]);
711 :
712 0 : fhPt[i] ->Add(entry->fhPt[i]);
713 0 : fhEta[i] ->Add(entry->fhEta[i]);
714 0 : fhTOFdistance[i] ->Add(entry->fhTOFdistance[i]);
715 : }
716 :
717 0 : fhCutStatistics ->Add(entry->fhCutStatistics);
718 0 : fhCutCorrelation ->Add(entry->fhCutCorrelation);
719 :
720 0 : count++;
721 0 : }
722 0 : return count+1;
723 0 : }
724 :
725 : void AliESDtrackCuts::SetMinNClustersTPCPtDep(TFormula *f1, Float_t ptmax)
726 : {
727 : //
728 : // Sets the pT dependent NClustersTPC cut
729 : //
730 :
731 0 : if(f1){
732 0 : delete f1CutMinNClustersTPCPtDep;
733 0 : f1CutMinNClustersTPCPtDep = (TFormula*)f1->Clone("f1CutMinNClustersTPCPtDep");
734 0 : }
735 0 : fCutMaxPtDepNClustersTPC=ptmax;
736 0 : }
737 :
738 : void AliESDtrackCuts::SetCutGeoNcrNcl(Float_t deadZoneWidth,Float_t cutGeoNcrNclLength, Float_t cutGeoNcrNclGeom1Pt, Float_t cutGeoNcrNclFractionNcr, Float_t cutGeoNcrNclFractionNcl){
739 : //
740 : // Set combination of the TPC geometrical cut, cut on number of crossed rows and cut on the clusters
741 : // The goal of cut is to avoid of the dead zone regions
742 : // - pointer to WIKI to put here
743 : // - JIRA discussion and links to presenttions - https://alice.its.cern.ch/jira/browse/ATO-233
744 : // Cut to be used for measurement where the homogenous coverage are not required
745 : // In case good momentum resolution and description of the matching efficiency needed, also for analysis reuiring homogenous phi coverage - this cut to be condered
746 : // Note this is repalcemnt of single cut
747 : //
748 0 : fDeadZoneWidth=deadZoneWidth;
749 0 : fCutGeoNcrNclLength=cutGeoNcrNclLength;
750 0 : fCutGeoNcrNclGeom1Pt=cutGeoNcrNclGeom1Pt;
751 0 : fCutGeoNcrNclFractionNcr=cutGeoNcrNclFractionNcr;
752 0 : fCutGeoNcrNclFractionNcl=cutGeoNcrNclFractionNcl;
753 0 : }
754 :
755 :
756 :
757 :
758 : //____________________________________________________________________
759 : AliESDtrackCuts* AliESDtrackCuts::GetStandardTPCOnlyTrackCuts()
760 : {
761 : // creates an AliESDtrackCuts object and fills it with standard (pre data-taking) values for TPC-only cuts
762 :
763 8 : AliInfoClass("Creating track cuts for TPC-only.");
764 :
765 4 : AliESDtrackCuts* esdTrackCuts = new AliESDtrackCuts;
766 :
767 4 : esdTrackCuts->SetMinNClustersTPC(50);
768 4 : esdTrackCuts->SetMaxChi2PerClusterTPC(4);
769 4 : esdTrackCuts->SetAcceptKinkDaughters(kFALSE);
770 :
771 4 : esdTrackCuts->SetMaxDCAToVertexZ(3.2);
772 4 : esdTrackCuts->SetMaxDCAToVertexXY(2.4);
773 4 : esdTrackCuts->SetDCAToVertex2D(kTRUE);
774 :
775 4 : return esdTrackCuts;
776 0 : }
777 :
778 : //____________________________________________________________________
779 : AliESDtrackCuts* AliESDtrackCuts::GetStandardITSTPCTrackCuts2009(Bool_t selPrimaries)
780 : {
781 : // creates an AliESDtrackCuts object and fills it with standard values for ITS-TPC cuts for pp 2009 data
782 :
783 0 : AliInfoClass("Creating track cuts for ITS+TPC (2009 definition).");
784 :
785 0 : AliESDtrackCuts* esdTrackCuts = new AliESDtrackCuts;
786 :
787 : // TPC
788 0 : esdTrackCuts->SetRequireTPCStandAlone(kTRUE); // to get chi2 and ncls of kTPCin
789 0 : esdTrackCuts->SetMinNClustersTPC(70);
790 0 : esdTrackCuts->SetMaxChi2PerClusterTPC(4);
791 0 : esdTrackCuts->SetAcceptKinkDaughters(kFALSE);
792 0 : esdTrackCuts->SetRequireTPCRefit(kTRUE);
793 : // ITS
794 0 : esdTrackCuts->SetRequireITSRefit(kTRUE);
795 0 : esdTrackCuts->SetClusterRequirementITS(AliESDtrackCuts::kSPD,
796 : AliESDtrackCuts::kAny);
797 0 : if(selPrimaries) {
798 : // 7*(0.0050+0.0060/pt^0.9)
799 0 : esdTrackCuts->SetMaxDCAToVertexXYPtDep("0.0350+0.0420/pt^0.9");
800 0 : esdTrackCuts->SetMaxChi2TPCConstrainedGlobal(36);
801 0 : }
802 0 : esdTrackCuts->SetMaxDCAToVertexZ(1.e6);
803 0 : esdTrackCuts->SetDCAToVertex2D(kFALSE);
804 0 : esdTrackCuts->SetRequireSigmaToVertex(kFALSE);
805 : //esdTrackCuts->SetEtaRange(-0.8,+0.8);
806 :
807 0 : esdTrackCuts->SetMaxChi2PerClusterITS(36);
808 :
809 0 : return esdTrackCuts;
810 0 : }
811 :
812 : //____________________________________________________________________
813 : AliESDtrackCuts* AliESDtrackCuts::GetStandardITSTPCTrackCuts2011(Bool_t selPrimaries, Int_t clusterCut)
814 : {
815 : // creates an AliESDtrackCuts object and fills it with standard values for ITS-TPC cuts for pp 2011 data
816 : // if clusterCut = 1, the cut on the number of clusters is replaced by
817 : // a cut on the number of crossed rows and on the ration crossed
818 : // rows/findable clusters
819 :
820 0 : AliInfoClass("Creating track cuts for ITS+TPC (2011 definition).");
821 :
822 0 : AliESDtrackCuts* esdTrackCuts = new AliESDtrackCuts;
823 :
824 : // TPC
825 0 : if(clusterCut == 0) esdTrackCuts->SetMinNClustersTPC(50);
826 0 : else if (clusterCut == 1) {
827 0 : esdTrackCuts->SetMinNCrossedRowsTPC(70);
828 0 : esdTrackCuts->SetMinRatioCrossedRowsOverFindableClustersTPC(0.8);
829 0 : }
830 : else {
831 0 : AliWarningClass(Form("Wrong value of the clusterCut parameter (%d), using cut on Nclusters",clusterCut));
832 0 : esdTrackCuts->SetMinNClustersTPC(50);
833 : }
834 0 : esdTrackCuts->SetMaxChi2PerClusterTPC(4);
835 0 : esdTrackCuts->SetAcceptKinkDaughters(kFALSE);
836 0 : esdTrackCuts->SetRequireTPCRefit(kTRUE);
837 : // ITS
838 0 : esdTrackCuts->SetRequireITSRefit(kTRUE);
839 0 : esdTrackCuts->SetClusterRequirementITS(AliESDtrackCuts::kSPD,
840 : AliESDtrackCuts::kAny);
841 0 : if(selPrimaries) {
842 : // 7*(0.0015+0.0050/pt^1.1)
843 0 : esdTrackCuts->SetMaxDCAToVertexXYPtDep("0.0105+0.0350/pt^1.1");
844 0 : esdTrackCuts->SetMaxChi2TPCConstrainedGlobal(36);
845 0 : }
846 0 : esdTrackCuts->SetMaxDCAToVertexZ(2);
847 0 : esdTrackCuts->SetDCAToVertex2D(kFALSE);
848 0 : esdTrackCuts->SetRequireSigmaToVertex(kFALSE);
849 :
850 0 : esdTrackCuts->SetMaxChi2PerClusterITS(36);
851 :
852 0 : return esdTrackCuts;
853 0 : }
854 :
855 : //____________________________________________________________________
856 : AliESDtrackCuts* AliESDtrackCuts::GetStandardITSTPCTrackCuts2015PbPb(Bool_t selPrimaries/*=kTRUE*/, Int_t clusterCut/*=1*/, Bool_t cutAcceptanceEdges/*=kTRUE*/, Bool_t removeDistortedRegions/*=kFALSE*/) {
857 : //
858 : // creates an AliESDtrackCuts object and fills it with standard values for ITS-TPC cuts for PbPb 2015 data
859 : // if clusterCut = 1, the cut on the number of clusters is replaced by
860 : // a cut on the number of crossed rows and on the ration crossed
861 : // rows/findable clusters
862 : //
863 0 : AliInfoClass("Creating track cuts for ITS+TPC (2015 Pb-Pb run definition).");
864 0 : AliWarningClass("THIS TRACK-CUT SET HAS NOT YET BEEN VALIDATED AND IS NOT YET FROZEN. TO BE USED FOR TESTING PURPOSES ONLY.");
865 : //
866 0 : AliESDtrackCuts* esdTrackCuts = new AliESDtrackCuts;
867 :
868 : // TPC
869 0 : if(clusterCut == 0) esdTrackCuts->SetMinNClustersTPC(50);
870 0 : else if (clusterCut == 1) {
871 0 : esdTrackCuts->SetMinNCrossedRowsTPC(70);
872 0 : esdTrackCuts->SetMinRatioCrossedRowsOverFindableClustersTPC(0.8);
873 0 : }
874 : else {
875 0 : AliWarningClass(Form("Wrong value of the clusterCut parameter (%d), using cut on Nclusters",clusterCut));
876 0 : esdTrackCuts->SetMinNClustersTPC(50);
877 : }
878 : //
879 0 : if (cutAcceptanceEdges) esdTrackCuts->SetCutGeoNcrNcl(2., 130., 1.5, 0.0, 0.0); // only dead zone and not clusters per length
880 0 : if (removeDistortedRegions) esdTrackCuts->SetCutOutDistortedRegionsTPC(kTRUE);
881 : //
882 0 : esdTrackCuts->SetMaxChi2PerClusterTPC(4);
883 0 : esdTrackCuts->SetAcceptKinkDaughters(kFALSE);
884 0 : esdTrackCuts->SetRequireTPCRefit(kTRUE);
885 : // ITS
886 0 : esdTrackCuts->SetRequireITSRefit(kTRUE);
887 0 : esdTrackCuts->SetClusterRequirementITS(AliESDtrackCuts::kSPD,
888 : AliESDtrackCuts::kAny);
889 0 : if(selPrimaries) {
890 : // 7*(0.0015+0.0050/pt^1.1)
891 0 : esdTrackCuts->SetMaxDCAToVertexXYPtDep("0.0105+0.0350/pt^1.1");
892 0 : esdTrackCuts->SetMaxChi2TPCConstrainedGlobal(36);
893 0 : }
894 0 : esdTrackCuts->SetMaxDCAToVertexZ(2);
895 0 : esdTrackCuts->SetDCAToVertex2D(kFALSE);
896 0 : esdTrackCuts->SetRequireSigmaToVertex(kFALSE);
897 :
898 0 : esdTrackCuts->SetMaxChi2PerClusterITS(36);
899 :
900 0 : return esdTrackCuts;
901 :
902 :
903 0 : }
904 :
905 :
906 : //____________________________________________________________________
907 : AliESDtrackCuts* AliESDtrackCuts::GetStandardITSTPCTrackCuts2010(Bool_t selPrimaries,Int_t clusterCut)
908 : {
909 : // creates an AliESDtrackCuts object and fills it with standard values for ITS-TPC cuts for pp 2010 data
910 : // if clusterCut = 1, the cut on the number of clusters is replaced by
911 : // a cut on the number of crossed rows and on the ration crossed
912 : // rows/findable clusters
913 :
914 0 : AliInfoClass("Creating track cuts for ITS+TPC (2010 definition).");
915 :
916 0 : AliESDtrackCuts* esdTrackCuts = new AliESDtrackCuts;
917 :
918 : // TPC
919 0 : if(clusterCut == 0) esdTrackCuts->SetMinNClustersTPC(70);
920 0 : else if (clusterCut == 1) {
921 0 : esdTrackCuts->SetMinNCrossedRowsTPC(70);
922 0 : esdTrackCuts->SetMinRatioCrossedRowsOverFindableClustersTPC(0.8);
923 0 : }
924 : else {
925 0 : AliWarningClass(Form("Wrong value of the clusterCut parameter (%d), using cut on Nclusters",clusterCut));
926 0 : esdTrackCuts->SetMinNClustersTPC(70);
927 : }
928 0 : esdTrackCuts->SetMaxChi2PerClusterTPC(4);
929 0 : esdTrackCuts->SetAcceptKinkDaughters(kFALSE);
930 0 : esdTrackCuts->SetRequireTPCRefit(kTRUE);
931 : // ITS
932 0 : esdTrackCuts->SetRequireITSRefit(kTRUE);
933 0 : esdTrackCuts->SetClusterRequirementITS(AliESDtrackCuts::kSPD,
934 : AliESDtrackCuts::kAny);
935 0 : if(selPrimaries) {
936 : // 7*(0.0026+0.0050/pt^1.01)
937 0 : esdTrackCuts->SetMaxDCAToVertexXYPtDep("0.0182+0.0350/pt^1.01");
938 0 : esdTrackCuts->SetMaxChi2TPCConstrainedGlobal(36);
939 0 : }
940 0 : esdTrackCuts->SetMaxDCAToVertexZ(2);
941 0 : esdTrackCuts->SetDCAToVertex2D(kFALSE);
942 0 : esdTrackCuts->SetRequireSigmaToVertex(kFALSE);
943 :
944 0 : esdTrackCuts->SetMaxChi2PerClusterITS(36);
945 :
946 0 : return esdTrackCuts;
947 0 : }
948 :
949 : //____________________________________________________________________
950 : AliESDtrackCuts* AliESDtrackCuts::GetStandardITSPureSATrackCuts2009(Bool_t selPrimaries, Bool_t useForPid)
951 : {
952 : // creates an AliESDtrackCuts object and fills it with standard values for ITS pure SA tracks
953 :
954 0 : AliESDtrackCuts* esdTrackCuts = new AliESDtrackCuts;
955 0 : esdTrackCuts->SetRequireITSPureStandAlone(kTRUE);
956 0 : esdTrackCuts->SetRequireITSRefit(kTRUE);
957 0 : esdTrackCuts->SetMinNClustersITS(4);
958 0 : esdTrackCuts->SetClusterRequirementITS(AliESDtrackCuts::kSPD,
959 : AliESDtrackCuts::kAny);
960 0 : esdTrackCuts->SetMaxChi2PerClusterITS(1.);
961 :
962 0 : if(selPrimaries) {
963 : // 7*(0.0085+0.0026/pt^1.55)
964 0 : esdTrackCuts->SetMaxDCAToVertexXYPtDep("0.0595+0.0182/pt^1.55");
965 0 : }
966 0 : if(useForPid){
967 0 : esdTrackCuts->SetRequireITSPid(kTRUE);
968 0 : }
969 0 : return esdTrackCuts;
970 0 : }
971 :
972 : //____________________________________________________________________
973 : AliESDtrackCuts* AliESDtrackCuts::GetStandardITSPureSATrackCuts2010(Bool_t selPrimaries, Bool_t useForPid)
974 : {
975 : // creates an AliESDtrackCuts object and fills it with standard values for ITS pure SA tracks - pp 2010
976 :
977 0 : AliESDtrackCuts* esdTrackCuts = new AliESDtrackCuts;
978 0 : esdTrackCuts->SetRequireITSPureStandAlone(kTRUE);
979 0 : esdTrackCuts->SetRequireITSRefit(kTRUE);
980 0 : esdTrackCuts->SetMinNClustersITS(4);
981 0 : esdTrackCuts->SetClusterRequirementITS(AliESDtrackCuts::kSPD,
982 : AliESDtrackCuts::kAny);
983 0 : esdTrackCuts->SetMaxChi2PerClusterITS(2.5);
984 :
985 0 : if(selPrimaries) {
986 : // 7*(0.0033+0.0045/pt^1.3)
987 0 : esdTrackCuts->SetMaxDCAToVertexXYPtDep("0.0231+0.0315/pt^1.3");
988 0 : }
989 0 : if(useForPid){
990 0 : esdTrackCuts->SetRequireITSPid(kTRUE);
991 0 : }
992 0 : return esdTrackCuts;
993 0 : }
994 :
995 : //____________________________________________________________________
996 : AliESDtrackCuts* AliESDtrackCuts::GetStandardITSSATrackCuts2009(Bool_t selPrimaries, Bool_t useForPid)
997 : {
998 : // creates an AliESDtrackCuts object and fills it with standard values for ITS pure SA tracks
999 :
1000 0 : AliESDtrackCuts* esdTrackCuts = new AliESDtrackCuts;
1001 0 : esdTrackCuts->SetRequireITSStandAlone(kTRUE);
1002 0 : esdTrackCuts->SetRequireITSPureStandAlone(kFALSE);
1003 0 : esdTrackCuts->SetRequireITSRefit(kTRUE);
1004 0 : esdTrackCuts->SetMinNClustersITS(4);
1005 0 : esdTrackCuts->SetClusterRequirementITS(AliESDtrackCuts::kSPD,
1006 : AliESDtrackCuts::kAny);
1007 0 : esdTrackCuts->SetMaxChi2PerClusterITS(1.);
1008 :
1009 0 : if(selPrimaries) {
1010 : // 7*(0.0085+0.0026/pt^1.55)
1011 0 : esdTrackCuts->SetMaxDCAToVertexXYPtDep("0.0595+0.0182/pt^1.55");
1012 0 : }
1013 0 : if(useForPid){
1014 0 : esdTrackCuts->SetRequireITSPid(kTRUE);
1015 0 : }
1016 0 : return esdTrackCuts;
1017 0 : }
1018 :
1019 : //____________________________________________________________________
1020 : AliESDtrackCuts* AliESDtrackCuts::GetStandardITSSATrackCuts2010(Bool_t selPrimaries, Bool_t useForPid)
1021 : {
1022 : // creates an AliESDtrackCuts object and fills it with standard values for ITS pure SA tracks --pp 2010
1023 :
1024 0 : AliESDtrackCuts* esdTrackCuts = new AliESDtrackCuts;
1025 0 : esdTrackCuts->SetRequireITSStandAlone(kTRUE);
1026 0 : esdTrackCuts->SetRequireITSPureStandAlone(kFALSE);
1027 0 : esdTrackCuts->SetRequireITSRefit(kTRUE);
1028 0 : esdTrackCuts->SetMinNClustersITS(4);
1029 0 : esdTrackCuts->SetClusterRequirementITS(AliESDtrackCuts::kSPD,
1030 : AliESDtrackCuts::kAny);
1031 0 : esdTrackCuts->SetMaxChi2PerClusterITS(2.5);
1032 :
1033 0 : if(selPrimaries) {
1034 : // 7*(0.0033+0.0045/pt^1.3)
1035 0 : esdTrackCuts->SetMaxDCAToVertexXYPtDep("0.0231+0.0315/pt^1.3");
1036 0 : }
1037 0 : if(useForPid){
1038 0 : esdTrackCuts->SetRequireITSPid(kTRUE);
1039 0 : }
1040 0 : return esdTrackCuts;
1041 0 : }
1042 :
1043 : //____________________________________________________________________
1044 : AliESDtrackCuts* AliESDtrackCuts::GetStandardITSSATrackCutsPbPb2010(Bool_t selPrimaries, Bool_t useForPid)
1045 : {
1046 : // creates an AliESDtrackCuts object and fills it with standard values for ITS pure SA tracks -- PbPb 2010
1047 :
1048 0 : AliESDtrackCuts* esdTrackCuts = GetStandardITSSATrackCuts2010(selPrimaries, useForPid);
1049 0 : esdTrackCuts->SetMaxNOfMissingITSPoints(1);
1050 :
1051 0 : return esdTrackCuts;
1052 : }
1053 : //____________________________________________________________________
1054 :
1055 : AliESDtrackCuts* AliESDtrackCuts::GetStandardV0DaughterCuts()
1056 : {
1057 : // creates a AliESDtrackCuts object and fills it with standard cuts for V0 daughters
1058 0 : AliESDtrackCuts* esdTrackCuts = new AliESDtrackCuts;
1059 0 : esdTrackCuts->SetRequireTPCRefit(kTRUE);
1060 0 : esdTrackCuts->SetMinNClustersTPC(70);
1061 0 : esdTrackCuts->SetAcceptKinkDaughters(kFALSE);
1062 0 : return esdTrackCuts;
1063 0 : }
1064 :
1065 : //____________________________________________________________________
1066 : Bool_t AliESDtrackCuts::IsTrackInDistortedTpcRegion(const AliESDtrack * esdTrack) {
1067 : //
1068 : // Returns kTRUE if the track crosses on of the regions in the TPC
1069 : // in which the field is distorted. For the time being, this is a
1070 : // purely geometric cut which will be further refined in the future
1071 : // based on NDregression analysis.
1072 : //
1073 : // For the time being, the cuts are hardwired, but will be probably
1074 : // moved to the OADB.
1075 : //
1076 0 : if (!esdTrack->GetInnerParam()) return kFALSE; // non-tpc tracks are not affected
1077 : //
1078 0 : Double_t eta = esdTrack->Eta();
1079 0 : Double_t phiIn = esdTrack->GetInnerParam()->Phi();
1080 : //
1081 0 : if (eta > 0) { // A-side
1082 0 : if (1.32 < phiIn && phiIn < 1.48) return kTRUE;
1083 0 : if (2.00 < phiIn && phiIn < 2.15) return kTRUE;
1084 0 : if (3.07 < phiIn && phiIn < 3.21) return kTRUE;
1085 : } else { // C-side
1086 0 : if (0.40 < phiIn && phiIn < 0.86) return kTRUE;
1087 0 : if (2.10 < phiIn && phiIn < 2.34) return kTRUE;
1088 0 : if (3.90 < phiIn && phiIn < 4.32) return kTRUE;
1089 : }
1090 : //
1091 0 : return kFALSE;
1092 0 : }
1093 :
1094 :
1095 : //____________________________________________________________________
1096 : Int_t AliESDtrackCuts::GetReferenceMultiplicity(const AliESDEvent* esd, Bool_t tpcOnly)
1097 : {
1098 : // Gets reference multiplicity following the standard cuts and a defined fiducial volume
1099 : // tpcOnly = kTRUE -> consider TPC-only tracks
1100 : // = kFALSE -> consider global tracks
1101 : //
1102 : // DEPRECATED Use GetReferenceMultiplicity with the enum as second argument instead
1103 :
1104 8 : if (!tpcOnly)
1105 : {
1106 0 : AliErrorClass("Not implemented for global tracks!");
1107 0 : return -1;
1108 : }
1109 :
1110 : static AliESDtrackCuts* esdTrackCuts = 0;
1111 8 : if (!esdTrackCuts)
1112 : {
1113 2 : esdTrackCuts = GetStandardTPCOnlyTrackCuts();
1114 2 : esdTrackCuts->SetEtaRange(-0.8, 0.8);
1115 2 : esdTrackCuts->SetPtRange(0.15);
1116 2 : }
1117 :
1118 8 : Int_t nTracks = esdTrackCuts->CountAcceptedTracks(esd);
1119 :
1120 : return nTracks;
1121 8 : }
1122 :
1123 : //____________________________________________________________________
1124 : Float_t AliESDtrackCuts::GetSigmaToVertex(const AliESDtrack* const esdTrack)
1125 : {
1126 : // Calculates the number of sigma to the vertex.
1127 :
1128 1316 : Float_t b[2];
1129 : Float_t bRes[2];
1130 658 : Float_t bCov[3];
1131 658 : esdTrack->GetImpactParameters(b,bCov);
1132 :
1133 1316 : if (bCov[0]<=0 || bCov[2]<=0) {
1134 0 : AliDebugClass(1, "Estimated b resolution lower or equal zero!");
1135 0 : bCov[0]=0; bCov[2]=0;
1136 0 : }
1137 658 : bRes[0] = TMath::Sqrt(bCov[0]);
1138 658 : bRes[1] = TMath::Sqrt(bCov[2]);
1139 :
1140 : // -----------------------------------
1141 : // How to get to a n-sigma cut?
1142 : //
1143 : // The accumulated statistics from 0 to d is
1144 : //
1145 : // -> Erf(d/Sqrt(2)) for a 1-dim gauss (d = n_sigma)
1146 : // -> 1 - Exp(-d**2) for a 2-dim gauss (d*d = dx*dx + dy*dy != n_sigma)
1147 : //
1148 : // It means that for a 2-dim gauss: n_sigma(d) = Sqrt(2)*ErfInv(1 - Exp((-d**2)/2)
1149 : // Can this be expressed in a different way?
1150 :
1151 1316 : if (bRes[0] == 0 || bRes[1] ==0)
1152 0 : return -1;
1153 :
1154 658 : Float_t d = TMath::Sqrt(TMath::Power(b[0]/bRes[0],2) + TMath::Power(b[1]/bRes[1],2));
1155 :
1156 : // work around precision problem
1157 : // if d is too big, TMath::Exp(...) gets 0, and TMath::ErfInverse(1) that should be infinite, gets 0 :(
1158 : // 1e-15 corresponds to nsigma ~ 7.7
1159 658 : if (TMath::Exp(-d * d / 2) < 1e-15)
1160 0 : return 1000;
1161 :
1162 658 : Float_t nSigma = TMath::ErfInverse(1 - TMath::Exp(-d * d / 2)) * TMath::Sqrt(2);
1163 : return nSigma;
1164 658 : }
1165 :
1166 :
1167 : //____________________________________________________________________
1168 : Float_t AliESDtrackCuts::GetSigmaToVertexVTrack(const AliVTrack* const vTrack)
1169 : {
1170 : // Calculates the number of sigma to the vertex.
1171 : // This method utilizes AliVTracks, so it can be used
1172 : // for ESDs or flat ESDs.
1173 :
1174 0 : Float_t b[2] = {0.,0.};
1175 : Float_t bRes[2] = {0.,0.};
1176 0 : Float_t bCov[3] = {0.,0.,0.};
1177 0 : vTrack->GetImpactParameters(b,bCov);
1178 :
1179 0 : if (bCov[0]<=0 || bCov[2]<=0) {
1180 0 : AliDebugClass(1, "Estimated b resolution lower or equal zero!");
1181 0 : bCov[0]=0; bCov[2]=0;
1182 0 : }
1183 0 : bRes[0] = TMath::Sqrt(bCov[0]);
1184 0 : bRes[1] = TMath::Sqrt(bCov[2]);
1185 :
1186 : // -----------------------------------
1187 : // How to get to a n-sigma cut?
1188 : //
1189 : // The accumulated statistics from 0 to d is
1190 : //
1191 : // -> Erf(d/Sqrt(2)) for a 1-dim gauss (d = n_sigma)
1192 : // -> 1 - Exp(-d**2) for a 2-dim gauss (d*d = dx*dx + dy*dy != n_sigma)
1193 : //
1194 : // It means that for a 2-dim gauss: n_sigma(d) = Sqrt(2)*ErfInv(1 - Exp((-d**2)/2)
1195 : // Can this be expressed in a different way?
1196 :
1197 0 : if (bRes[0] == 0 || bRes[1] ==0)
1198 0 : return -1;
1199 :
1200 0 : Float_t d = TMath::Sqrt(TMath::Power(b[0]/bRes[0],2) + TMath::Power(b[1]/bRes[1],2));
1201 :
1202 : // work around precision problem
1203 : // if d is too big, TMath::Exp(...) gets 0, and TMath::ErfInverse(1) that should be infinite, gets 0 :(
1204 : // 1e-15 corresponds to nsigma ~ 7.7
1205 0 : if (TMath::Exp(-d * d / 2) < 1e-15)
1206 0 : return 1000;
1207 :
1208 0 : Float_t nSigma = TMath::ErfInverse(1 - TMath::Exp(-d * d / 2)) * TMath::Sqrt(2);
1209 : return nSigma;
1210 0 : }
1211 :
1212 :
1213 : void AliESDtrackCuts::EnableNeededBranches(TTree* tree)
1214 : {
1215 : // enables the branches needed by AcceptTrack, for a list see comment of AcceptTrack
1216 :
1217 0 : tree->SetBranchStatus("fTracks.fFlags", 1);
1218 0 : tree->SetBranchStatus("fTracks.fITSncls", 1);
1219 0 : tree->SetBranchStatus("fTracks.fTPCncls", 1);
1220 0 : tree->SetBranchStatus("fTracks.fITSchi2", 1);
1221 0 : tree->SetBranchStatus("fTracks.fTPCchi2", 1);
1222 0 : tree->SetBranchStatus("fTracks.fC*", 1);
1223 0 : tree->SetBranchStatus("fTracks.fD", 1);
1224 0 : tree->SetBranchStatus("fTracks.fZ", 1);
1225 0 : tree->SetBranchStatus("fTracks.fCdd", 1);
1226 0 : tree->SetBranchStatus("fTracks.fCdz", 1);
1227 0 : tree->SetBranchStatus("fTracks.fCzz", 1);
1228 0 : tree->SetBranchStatus("fTracks.fP*", 1);
1229 0 : tree->SetBranchStatus("fTracks.fR*", 1);
1230 0 : tree->SetBranchStatus("fTracks.fKinkIndexes*", 1);
1231 0 : }
1232 :
1233 : //____________________________________________________________________
1234 : Bool_t AliESDtrackCuts::AcceptTrack(const AliESDtrack* esdTrack)
1235 : {
1236 : //
1237 : // figure out if the tracks survives all the track cuts defined
1238 : //
1239 : // the different quality parameter and kinematic values are first
1240 : // retrieved from the track. then it is found out what cuts the
1241 : // track did not survive and finally the cuts are imposed.
1242 :
1243 : // this function needs the following branches:
1244 : // fTracks.fFlags
1245 : // fTracks.fITSncls
1246 : // fTracks.fTPCncls
1247 : // fTracks.fITSchi2
1248 : // fTracks.fTPCchi2
1249 : // fTracks.fC //GetExternalCovariance
1250 : // fTracks.fD //GetImpactParameters
1251 : // fTracks.fZ //GetImpactParameters
1252 : // fTracks.fCdd //GetImpactParameters
1253 : // fTracks.fCdz //GetImpactParameters
1254 : // fTracks.fCzz //GetImpactParameters
1255 : // fTracks.fP //GetPxPyPz
1256 : // fTracks.fR //GetMass
1257 : // fTracks.fP //GetMass
1258 : // fTracks.fKinkIndexes
1259 : //
1260 : // esdEvent is only required for the MaxChi2TPCConstrainedVsGlobal
1261 :
1262 2088 : UInt_t status = esdTrack->GetStatus();
1263 :
1264 : // getting quality parameters from the ESD track
1265 1044 : Int_t nClustersITS = esdTrack->GetITSclusters(0);
1266 : Int_t nClustersTPC = -1;
1267 1044 : if(fCutRequireTPCStandAlone) {
1268 0 : nClustersTPC = esdTrack->GetTPCNclsIter1();
1269 0 : }
1270 : else {
1271 1044 : nClustersTPC = esdTrack->GetTPCclusters(0);
1272 : }
1273 :
1274 : //Pt dependent NClusters Cut
1275 1044 : if(f1CutMinNClustersTPCPtDep) {
1276 0 : if(esdTrack->Pt()<fCutMaxPtDepNClustersTPC)
1277 0 : fCutMinNClusterTPC = (Int_t)(f1CutMinNClustersTPCPtDep->Eval(esdTrack->Pt()));
1278 : else
1279 0 : fCutMinNClusterTPC = (Int_t)(f1CutMinNClustersTPCPtDep->Eval(fCutMaxPtDepNClustersTPC));
1280 : }
1281 :
1282 1044 : Float_t nCrossedRowsTPC = esdTrack->GetTPCCrossedRows();
1283 : Float_t ratioCrossedRowsOverFindableClustersTPC = 1.0;
1284 1044 : if (esdTrack->GetTPCNclsF()>0) {
1285 894 : ratioCrossedRowsOverFindableClustersTPC = nCrossedRowsTPC / esdTrack->GetTPCNclsF();
1286 894 : }
1287 :
1288 1044 : Int_t nClustersTPCShared = esdTrack->GetTPCnclsS();
1289 : Float_t fracClustersTPCShared = -1.;
1290 :
1291 : Float_t chi2PerClusterITS = -1;
1292 : Float_t chi2PerClusterTPC = -1;
1293 1044 : if (nClustersITS!=0)
1294 902 : chi2PerClusterITS = esdTrack->GetITSchi2()/Float_t(nClustersITS);
1295 1044 : if (nClustersTPC!=0) {
1296 894 : if(fCutRequireTPCStandAlone) {
1297 0 : chi2PerClusterTPC = esdTrack->GetTPCchi2Iter1()/Float_t(nClustersTPC);
1298 0 : } else {
1299 894 : chi2PerClusterTPC = esdTrack->GetTPCchi2()/Float_t(nClustersTPC);
1300 : }
1301 894 : fracClustersTPCShared = Float_t(nClustersTPCShared)/Float_t(nClustersTPC);
1302 894 : }
1303 :
1304 1044 : Double_t extCov[15];
1305 1044 : esdTrack->GetExternalCovariance(extCov);
1306 :
1307 1044 : Float_t b[2];
1308 1044 : Float_t bCov[3];
1309 1044 : esdTrack->GetImpactParameters(b,bCov);
1310 2088 : if (bCov[0]<=0 || bCov[2]<=0) {
1311 0 : AliDebug(1, "Estimated b resolution lower or equal zero!");
1312 0 : bCov[0]=0; bCov[2]=0;
1313 0 : }
1314 :
1315 :
1316 : // set pt-dependent DCA cuts, if requested
1317 1044 : SetPtDepDCACuts(esdTrack->Pt());
1318 :
1319 :
1320 1044 : Float_t dcaToVertexXY = b[0];
1321 1044 : Float_t dcaToVertexZ = b[1];
1322 :
1323 : Float_t dcaToVertex = -1;
1324 :
1325 2088 : if (fCutDCAToVertex2D)
1326 : {
1327 1343 : dcaToVertex = TMath::Sqrt(dcaToVertexXY*dcaToVertexXY/fCutMaxDCAToVertexXY/fCutMaxDCAToVertexXY + dcaToVertexZ*dcaToVertexZ/fCutMaxDCAToVertexZ/fCutMaxDCAToVertexZ);
1328 299 : }
1329 : else
1330 745 : dcaToVertex = TMath::Sqrt(dcaToVertexXY*dcaToVertexXY + dcaToVertexZ*dcaToVertexZ);
1331 :
1332 : // getting the kinematic variables of the track
1333 : // (assuming the mass is known)
1334 1044 : Double_t p[3];
1335 1044 : esdTrack->GetPxPyPz(p);
1336 :
1337 : // Changed from float to double to prevent rounding errors leading to negative
1338 : // log arguments (M.G.)
1339 1044 : Double_t momentum = TMath::Sqrt(p[0]*p[0] + p[1]*p[1] + p[2]*p[2]);
1340 1044 : Double_t pt = TMath::Sqrt(p[0]*p[0] + p[1]*p[1]);
1341 1044 : Double_t mass = esdTrack->GetMass();
1342 1044 : Double_t energy = TMath::Sqrt(mass*mass + momentum*momentum);
1343 :
1344 : //y-eta related calculations
1345 : Float_t eta = -100.;
1346 : Float_t y = -100.;
1347 1044 : if((momentum != TMath::Abs(p[2]))&&(momentum != 0))
1348 1044 : eta = 0.5*TMath::Log((momentum + p[2])/(momentum - p[2]));
1349 1044 : if((energy != TMath::Abs(p[2]))&&(energy != 0))
1350 1044 : y = 0.5*TMath::Log((energy + p[2])/(energy - p[2]));
1351 :
1352 1044 : if (extCov[14] < 0)
1353 : {
1354 0 : AliWarning(Form("GetSigma1Pt2() returns negative value for external covariance matrix element fC[14]: %f. Corrupted track information, track will not be accepted!", extCov[14]));
1355 0 : return kFALSE;
1356 : }
1357 1044 : Float_t relUncertainty1Pt = TMath::Sqrt(extCov[14])*pt;
1358 :
1359 : //########################################################################
1360 : // cut the track?
1361 :
1362 1044 : Bool_t cuts[kNCuts];
1363 96048 : for (Int_t i=0; i<kNCuts; i++) cuts[i]=kFALSE;
1364 :
1365 : // track quality cuts
1366 1447 : if (fCutRequireTPCRefit && (status&AliESDtrack::kTPCrefit)==0)
1367 16 : cuts[0]=kTRUE;
1368 1044 : if (fCutRequireTPCStandAlone && (status&AliESDtrack::kTPCin)==0)
1369 0 : cuts[1]=kTRUE;
1370 1338 : if (fCutRequireITSRefit && (status&AliESDtrack::kITSrefit)==0)
1371 18 : cuts[2]=kTRUE;
1372 1044 : if (nClustersTPC<fCutMinNClusterTPC)
1373 40 : cuts[3]=kTRUE;
1374 1044 : if (nClustersITS<fCutMinNClusterITS)
1375 0 : cuts[4]=kTRUE;
1376 1044 : if (chi2PerClusterTPC>fCutMaxChi2PerClusterTPC)
1377 0 : cuts[5]=kTRUE;
1378 1044 : if (chi2PerClusterITS>fCutMaxChi2PerClusterITS)
1379 0 : cuts[6]=kTRUE;
1380 1044 : if (extCov[0] > fCutMaxC11)
1381 22 : cuts[7]=kTRUE;
1382 1044 : if (extCov[2] > fCutMaxC22)
1383 18 : cuts[8]=kTRUE;
1384 1044 : if (extCov[5] > fCutMaxC33)
1385 0 : cuts[9]=kTRUE;
1386 1044 : if (extCov[9] > fCutMaxC44)
1387 0 : cuts[10]=kTRUE;
1388 1044 : if (extCov[14] > fCutMaxC55)
1389 2 : cuts[11]=kTRUE;
1390 :
1391 : // cut 12 and 13 see below
1392 :
1393 1589 : if (!fCutAcceptKinkDaughters && esdTrack->GetKinkIndex(0)>0)
1394 8 : cuts[14]=kTRUE;
1395 : // track kinematics cut
1396 2088 : if((momentum < fPMin) || (momentum > fPMax))
1397 0 : cuts[15]=kTRUE;
1398 2072 : if((pt < fPtMin) || (pt > fPtMax))
1399 16 : cuts[16] = kTRUE;
1400 2088 : if((p[0] < fPxMin) || (p[0] > fPxMax))
1401 0 : cuts[17] = kTRUE;
1402 2088 : if((p[1] < fPyMin) || (p[1] > fPyMax))
1403 0 : cuts[18] = kTRUE;
1404 2088 : if((p[2] < fPzMin) || (p[2] > fPzMax))
1405 0 : cuts[19] = kTRUE;
1406 2088 : if((eta < fEtaMin) || (eta > fEtaMax))
1407 0 : cuts[20] = kTRUE;
1408 2088 : if((y < fRapMin) || (y > fRapMax))
1409 0 : cuts[21] = kTRUE;
1410 1044 : if (fCutDCAToVertex2D && dcaToVertex > 1)
1411 59 : cuts[22] = kTRUE;
1412 1789 : if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexXY) > fCutMaxDCAToVertexXY)
1413 18 : cuts[23] = kTRUE;
1414 1789 : if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexZ) > fCutMaxDCAToVertexZ)
1415 0 : cuts[24] = kTRUE;
1416 1343 : if (fCutDCAToVertex2D && fCutMinDCAToVertexXY > 0 && fCutMinDCAToVertexZ > 0 && dcaToVertexXY*dcaToVertexXY/fCutMinDCAToVertexXY/fCutMinDCAToVertexXY + dcaToVertexZ*dcaToVertexZ/fCutMinDCAToVertexZ/fCutMinDCAToVertexZ < 1)
1417 0 : cuts[25] = kTRUE;
1418 1789 : if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexXY) < fCutMinDCAToVertexXY)
1419 0 : cuts[26] = kTRUE;
1420 1789 : if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexZ) < fCutMinDCAToVertexZ)
1421 0 : cuts[27] = kTRUE;
1422 :
1423 8352 : for (Int_t i = 0; i < 3; i++) {
1424 3132 : if(!(esdTrack->GetStatus()&AliESDtrack::kITSupg)) { // current ITS
1425 3132 : cuts[28+i] = !CheckITSClusterRequirement(fCutClusterRequirementITS[i], esdTrack->HasPointOnITSLayer(i*2), esdTrack->HasPointOnITSLayer(i*2+1));
1426 3132 : } else { // upgraded ITS (7 layers)
1427 : // at the moment, for L012 the layers 12 are considered together
1428 0 : if(i==0) { // L012
1429 0 : cuts[28+i] = !CheckITSClusterRequirement(fCutClusterRequirementITS[i], esdTrack->HasPointOnITSLayer(0), (esdTrack->HasPointOnITSLayer(1))&(esdTrack->HasPointOnITSLayer(2)));
1430 0 : } else { // L34 or L56
1431 0 : cuts[28+i] = !CheckITSClusterRequirement(fCutClusterRequirementITS[i], esdTrack->HasPointOnITSLayer(i*2+1), esdTrack->HasPointOnITSLayer(i*2+2));
1432 : }
1433 : }
1434 : }
1435 :
1436 1931 : if(fCutRequireITSStandAlone || fCutRequireITSpureSA){
1437 265 : if ((status & AliESDtrack::kITSin) == 0 || (status & AliESDtrack::kTPCin)){
1438 : // TPC tracks
1439 141 : cuts[31] = kTRUE;
1440 141 : }else{
1441 : // ITS standalone tracks
1442 32 : if(fCutRequireITSStandAlone && !fCutRequireITSpureSA){
1443 16 : if(status & AliESDtrack::kITSpureSA) cuts[31] = kTRUE;
1444 0 : }else if(fCutRequireITSpureSA){
1445 0 : if(!(status & AliESDtrack::kITSpureSA)) cuts[31] = kTRUE;
1446 : }
1447 : }
1448 : }
1449 :
1450 1044 : if (relUncertainty1Pt > fCutMaxRel1PtUncertainty)
1451 0 : cuts[32] = kTRUE;
1452 :
1453 1044 : if (!fCutAcceptSharedTPCClusters && nClustersTPCShared!=0)
1454 0 : cuts[33] = kTRUE;
1455 :
1456 1044 : if (fracClustersTPCShared > fCutMaxFractionSharedTPCClusters)
1457 0 : cuts[34] = kTRUE;
1458 :
1459 : Int_t nITSPointsForPid=0;
1460 1044 : UChar_t clumap=esdTrack->GetITSClusterMap();
1461 10440 : for(Int_t i=2; i<6; i++){
1462 7686 : if(clumap&(1<<i)) ++nITSPointsForPid;
1463 : }
1464 1044 : if(fCutRequireITSPid && nITSPointsForPid<3) cuts[35] = kTRUE;
1465 :
1466 :
1467 1044 : if (nCrossedRowsTPC<fCutMinNCrossedRowsTPC)
1468 0 : cuts[36]=kTRUE;
1469 1044 : if (ratioCrossedRowsOverFindableClustersTPC<fCutMinRatioCrossedRowsOverFindableClustersTPC)
1470 0 : cuts[37]=kTRUE;
1471 :
1472 : Int_t nMissITSpts=0;
1473 1044 : Int_t idet,statusLay;
1474 1044 : Float_t xloc,zloc;
1475 14616 : for(Int_t iLay=0; iLay<6; iLay++){
1476 6264 : Bool_t retc=esdTrack->GetITSModuleIndexInfo(iLay,idet,statusLay,xloc,zloc);
1477 6318 : if(retc && statusLay==5) ++nMissITSpts;
1478 : }
1479 1044 : if(nMissITSpts>fCutMaxMissingITSPoints) cuts[38] = kTRUE;
1480 :
1481 : //kTOFout
1482 1044 : if (fCutRequireTOFout && (status&AliESDtrack::kTOFout)==0)
1483 0 : cuts[40]=kTRUE;
1484 :
1485 : // TOF signal Dz cut
1486 : Float_t dxTOF = 999.; //esdTrack->GetTOFsignalDx();
1487 : Float_t dzTOF = 999.; //esdTrack->GetTOFsignalDz();
1488 1044 : if (fFlagCutTOFdistance && (esdTrack->GetStatus() & AliESDtrack::kTOFout) == AliESDtrack::kTOFout){ // applying the TOF distance cut only if requested, and only on tracks that reached the TOF and where associated with a TOF hit
1489 0 : dxTOF = esdTrack->GetTOFsignalDx();
1490 0 : dzTOF = esdTrack->GetTOFsignalDz();
1491 0 : if (fgBeamTypeFlag < 0) { // the check on the beam type was not done yet
1492 0 : const AliESDEvent* event = esdTrack->GetESDEvent();
1493 0 : if (event){
1494 0 : TString beamTypeESD = event->GetBeamType();
1495 0 : AliDebug(2,Form("Beam type from ESD event = %s",beamTypeESD.Data()));
1496 0 : if (beamTypeESD.CompareTo("A-A",TString::kIgnoreCase) == 0){ // we are in PbPb collisions --> fgBeamTypeFlag will be set to 1, to apply the cut on TOF signal Dz
1497 0 : fgBeamTypeFlag = 1;
1498 0 : }
1499 : else { // we are NOT in PbPb collisions --> fgBeamTypeFlag will be set to 0, to NOT apply the cu6 on TOF signal Dz
1500 0 : fgBeamTypeFlag = 0;
1501 : }
1502 0 : }
1503 : else{
1504 0 : AliFatal("Beam type not available, but it is needed to apply the TOF cut!");
1505 : }
1506 0 : }
1507 :
1508 0 : if (fgBeamTypeFlag == 1){ // we are in PbPb collisions --> apply the cut on TOF signal Dz
1509 0 : Float_t radiusTOF = TMath::Sqrt(dxTOF*dxTOF + dzTOF*dzTOF);
1510 0 : AliDebug(3,Form("TOF check (with fCutTOFdistance = %f) --> dx = %f, dz = %f, radius = %f", fCutTOFdistance, dxTOF, dzTOF, radiusTOF));
1511 0 : if (radiusTOF > fCutTOFdistance){
1512 0 : AliDebug(2, Form("************* the radius is outside the range! %f > %f, the track will be skipped", radiusTOF, fCutTOFdistance));
1513 0 : cuts[41] = kTRUE;
1514 0 : }
1515 0 : }
1516 : }
1517 :
1518 : Bool_t cut=kFALSE;
1519 96048 : for (Int_t i=0; i<kNCuts; i++)
1520 47356 : if (cuts[i]) {cut = kTRUE;}
1521 :
1522 : // for performance evaluate the CPU intensive cuts only when the others have passed, and when they are requested
1523 : Double_t chi2TPCConstrainedVsGlobal = -2;
1524 : Float_t nSigmaToVertex = -2;
1525 1044 : if (!cut)
1526 : {
1527 : // getting the track to vertex parameters
1528 768 : if (fCutSigmaToVertexRequired)
1529 : {
1530 658 : nSigmaToVertex = GetSigmaToVertex(esdTrack);
1531 658 : if (nSigmaToVertex > fCutNsigmaToVertex && fCutSigmaToVertexRequired)
1532 : {
1533 0 : cuts[12] = kTRUE;
1534 : cut = kTRUE;
1535 0 : }
1536 : // if n sigma could not be calculated
1537 658 : if (nSigmaToVertex<0 && fCutSigmaToVertexRequired)
1538 : {
1539 0 : cuts[13] = kTRUE;
1540 : cut = kTRUE;
1541 0 : }
1542 : }
1543 :
1544 : // max chi2 TPC constrained vs global track only if track passed the other cut
1545 768 : if (fCutMaxChi2TPCConstrainedVsGlobal < 1e9)
1546 : {
1547 0 : const AliESDEvent* esdEvent = esdTrack->GetESDEvent();
1548 :
1549 0 : if (!esdEvent)
1550 0 : AliFatal("fCutMaxChi2TPCConstrainedVsGlobal set but ESD event not set in AliESDTrack. Use AliESDTrack::SetESDEvent before calling AliESDtrackCuts.");
1551 :
1552 : // get vertex
1553 : const AliESDVertex* vertex = 0;
1554 0 : if (fCutMaxChi2TPCConstrainedVsGlobalVertexType & kVertexTracks)
1555 0 : vertex = esdEvent->GetPrimaryVertexTracks();
1556 :
1557 0 : if ((!vertex || !vertex->GetStatus()) && fCutMaxChi2TPCConstrainedVsGlobalVertexType & kVertexSPD)
1558 0 : vertex = esdEvent->GetPrimaryVertexSPD();
1559 :
1560 0 : if ((!vertex || !vertex->GetStatus()) && fCutMaxChi2TPCConstrainedVsGlobalVertexType & kVertexTPC)
1561 0 : vertex = esdEvent->GetPrimaryVertexTPC();
1562 :
1563 0 : if (vertex->GetStatus())
1564 0 : chi2TPCConstrainedVsGlobal = esdTrack->GetChi2TPCConstrainedVsGlobal(vertex);
1565 :
1566 0 : if (chi2TPCConstrainedVsGlobal < 0 || chi2TPCConstrainedVsGlobal > fCutMaxChi2TPCConstrainedVsGlobal)
1567 : {
1568 0 : cuts[39] = kTRUE;
1569 : cut = kTRUE;
1570 0 : }
1571 0 : }
1572 :
1573 : // max length in active volume
1574 : Float_t lengthInActiveZoneTPC = -1;
1575 1536 : if (fCutMinLengthActiveVolumeTPC > 1. || fCutGeoNcrNclLength>0) { // do the calculation only if needed to save cpu-time
1576 0 : if (esdTrack->GetESDEvent()) {
1577 : //
1578 : const Double_t kMaxZ=220;
1579 0 : if (esdTrack->GetInnerParam()) lengthInActiveZoneTPC = esdTrack->GetLengthInActiveZone(1, fDeadZoneWidth, kMaxZ, esdTrack->GetESDEvent()->GetMagneticField());
1580 : //
1581 0 : if (lengthInActiveZoneTPC < fCutMinLengthActiveVolumeTPC ) {
1582 0 : cuts[42] = kTRUE;
1583 : cut = kTRUE;
1584 0 : }
1585 0 : if (fCutGeoNcrNclLength>0){
1586 0 : Double_t cutGeoNcrNclLength=fCutGeoNcrNclLength-TMath::Power(TMath::Abs(esdTrack->GetSigned1Pt()),fCutGeoNcrNclGeom1Pt);
1587 : Bool_t isOK=kTRUE;
1588 0 : if (lengthInActiveZoneTPC<cutGeoNcrNclLength) isOK=kFALSE;
1589 0 : if (nCrossedRowsTPC<fCutGeoNcrNclFractionNcr*cutGeoNcrNclLength) isOK=kFALSE;
1590 0 : if (esdTrack->GetTPCncls()<fCutGeoNcrNclFractionNcl*cutGeoNcrNclLength) isOK=kFALSE;
1591 :
1592 0 : if(!isOK) {
1593 0 : cuts[43] = kTRUE;
1594 : cut = kTRUE;
1595 0 : }
1596 0 : }
1597 0 : }
1598 : }
1599 :
1600 : // track in distorted TPC region
1601 768 : if (fCutOutDistortedRegionTPC) {
1602 0 : if (IsTrackInDistortedTpcRegion(esdTrack)) {
1603 0 : cuts[44] = kTRUE;
1604 : cut = kTRUE;
1605 0 : }
1606 : }
1607 768 : }
1608 :
1609 : //########################################################################
1610 : // filling histograms
1611 1044 : if (fHistogramsOn) {
1612 0 : fhCutStatistics->Fill(fhCutStatistics->GetBinCenter(fhCutStatistics->GetXaxis()->FindBin("n tracks")));
1613 0 : if (cut)
1614 0 : fhCutStatistics->Fill(fhCutStatistics->GetBinCenter(fhCutStatistics->GetXaxis()->FindBin("n cut tracks")));
1615 :
1616 0 : for (Int_t i=0; i<kNCuts; i++) {
1617 0 : if (fhCutStatistics->GetXaxis()->FindBin(fgkCutNames[i]) < 1)
1618 0 : AliFatal(Form("Inconsistency! Cut %d with name %s not found", i, fgkCutNames[i]));
1619 :
1620 0 : if (cuts[i])
1621 0 : fhCutStatistics->Fill(fhCutStatistics->GetBinCenter(fhCutStatistics->GetXaxis()->FindBin(fgkCutNames[i])));
1622 :
1623 0 : for (Int_t j=i; j<kNCuts; j++) {
1624 0 : if (cuts[i] && cuts[j]) {
1625 0 : Float_t xC = fhCutCorrelation->GetXaxis()->GetBinCenter(fhCutCorrelation->GetXaxis()->FindBin(fgkCutNames[i]));
1626 0 : Float_t yC = fhCutCorrelation->GetYaxis()->GetBinCenter(fhCutCorrelation->GetYaxis()->FindBin(fgkCutNames[j]));
1627 0 : fhCutCorrelation->Fill(xC, yC);
1628 0 : }
1629 : }
1630 : }
1631 0 : }
1632 :
1633 : // now we loop over the filling of the histograms twice: once "before" the cut, once "after"
1634 : // the code is not in a function due to too many local variables that would need to be passed
1635 :
1636 5928 : for (Int_t id = 0; id < 2; id++)
1637 : {
1638 : // id = 0 --> before cut
1639 : // id = 1 --> after cut
1640 :
1641 1812 : if (fHistogramsOn)
1642 : {
1643 0 : fhNClustersITS[id]->Fill(nClustersITS);
1644 0 : fhNClustersTPC[id]->Fill(nClustersTPC);
1645 0 : fhNSharedClustersTPC[id]->Fill(nClustersTPCShared);
1646 0 : fhNCrossedRowsTPC[id]->Fill(nCrossedRowsTPC);
1647 0 : fhRatioCrossedRowsOverFindableClustersTPC[id]->Fill(ratioCrossedRowsOverFindableClustersTPC);
1648 0 : fhChi2PerClusterITS[id]->Fill(chi2PerClusterITS);
1649 0 : fhChi2PerClusterTPC[id]->Fill(chi2PerClusterTPC);
1650 0 : fhChi2TPCConstrainedVsGlobal[id]->Fill(chi2TPCConstrainedVsGlobal);
1651 0 : fhNClustersForITSPID[id]->Fill(nITSPointsForPid);
1652 0 : fhNMissingITSPoints[id]->Fill(nMissITSpts);
1653 :
1654 0 : fhC11[id]->Fill(extCov[0]);
1655 0 : fhC22[id]->Fill(extCov[2]);
1656 0 : fhC33[id]->Fill(extCov[5]);
1657 0 : fhC44[id]->Fill(extCov[9]);
1658 0 : fhC55[id]->Fill(extCov[14]);
1659 :
1660 0 : fhRel1PtUncertainty[id]->Fill(relUncertainty1Pt);
1661 :
1662 0 : fhPt[id]->Fill(pt);
1663 0 : fhEta[id]->Fill(eta);
1664 0 : fhTOFdistance[id]->Fill(dxTOF, dzTOF);
1665 :
1666 : Float_t bRes[2];
1667 0 : bRes[0] = TMath::Sqrt(bCov[0]);
1668 0 : bRes[1] = TMath::Sqrt(bCov[2]);
1669 :
1670 0 : fhDZ[id]->Fill(b[1]);
1671 0 : fhDXY[id]->Fill(b[0]);
1672 0 : fhDXYDZ[id]->Fill(dcaToVertex);
1673 0 : fhDXYvsDZ[id]->Fill(b[1],b[0]);
1674 :
1675 0 : if (bRes[0]!=0 && bRes[1]!=0) {
1676 0 : fhDZNormalized[id]->Fill(b[1]/bRes[1]);
1677 0 : fhDXYNormalized[id]->Fill(b[0]/bRes[0]);
1678 0 : fhDXYvsDZNormalized[id]->Fill(b[1]/bRes[1], b[0]/bRes[0]);
1679 0 : fhNSigmaToVertex[id]->Fill(nSigmaToVertex);
1680 0 : }
1681 0 : }
1682 :
1683 : // cut the track
1684 1812 : if (cut)
1685 276 : return kFALSE;
1686 : }
1687 :
1688 768 : return kTRUE;
1689 2088 : }
1690 :
1691 : //____________________________________________________________________
1692 : Bool_t AliESDtrackCuts::AcceptVTrack(const AliVTrack* vTrack)
1693 : {
1694 : //
1695 : // figure out if the tracks survives all the track cuts defined
1696 : // Wrapper function designed to handle either AliESDtrack or
1697 : // AliFlatESDTrack objects using the AliVTrack interface as an
1698 : // intermediary. AliESDtracks will be shunted to AcceptTrack,
1699 : // while AliFlatESDtracks will continue in AcceptVTrack using
1700 : // a smaller subset of cuts (because many of the cuts can't be
1701 : // tested using the limited members contained within
1702 : // AliFlatESDtrack)
1703 : //
1704 : // the different quality parameter and kinematic values are first
1705 : // retrieved from the track. then it is found out what cuts the
1706 : // track did not survive and finally the cuts are imposed.
1707 :
1708 : // this function needs the following branches:
1709 : // fTracks.fFlags
1710 : // fTracks.fITSncls
1711 : // fTracks.fTPCncls
1712 : // fTracks.fITSchi2
1713 : // fTracks.fTPCchi2
1714 : // fTracks.fC //GetExternalCovariance
1715 : // fTracks.fD //GetImpactParameters
1716 : // fTracks.fZ //GetImpactParameters
1717 : // fTracks.fCdd //GetImpactParameters
1718 : // fTracks.fCdz //GetImpactParameters
1719 : // fTracks.fCzz //GetImpactParameters
1720 : // fTracks.fP //GetPxPyPz
1721 : // fTracks.fR //GetMass
1722 : // fTracks.fP //GetMass
1723 : // fTracks.fKinkIndexes
1724 : //
1725 : // esdEvent is only required for the MaxChi2TPCConstrainedVsGlobal
1726 :
1727 : //Check if the track is an AliESDtrack. If so, pass it to
1728 : //AcceptTrack() to use the full set of cuts
1729 0 : if(const AliESDtrack *esdTrack = dynamic_cast<const AliESDtrack*>(vTrack)){
1730 0 : return AcceptTrack(esdTrack);
1731 : }
1732 :
1733 : //The track is not an AliESDtrack. Perform a more limited
1734 : //set of cuts
1735 :
1736 0 : UInt_t status = vTrack->GetStatus();
1737 :
1738 : // getting quality parameters from the ESD track
1739 : // Int_t nClustersITS = vTrack->GetITSclusters(0);
1740 0 : Int_t nClustersITS = vTrack->GetNumberOfITSClusters();
1741 : // Int_t nClustersTPC = -1;
1742 : // if(fCutRequireTPCStandAlone) {
1743 : // nClustersTPC = vTrack->GetTPCNclsIter1();
1744 : // }
1745 : // else {
1746 : // nClustersTPC = vTrack->GetTPCclusters(0);
1747 : // }
1748 0 : Int_t nClustersTPC = vTrack->GetTPCNcls();
1749 :
1750 : //Pt dependent NClusters Cut
1751 0 : if(f1CutMinNClustersTPCPtDep) {
1752 0 : if(vTrack->Pt()<fCutMaxPtDepNClustersTPC)
1753 0 : fCutMinNClusterTPC = (Int_t)(f1CutMinNClustersTPCPtDep->Eval(vTrack->Pt()));
1754 : else
1755 0 : fCutMinNClusterTPC = (Int_t)(f1CutMinNClustersTPCPtDep->Eval(fCutMaxPtDepNClustersTPC));
1756 : }
1757 :
1758 0 : Float_t nCrossedRowsTPC = vTrack->GetTPCCrossedRows();
1759 : Float_t ratioCrossedRowsOverFindableClustersTPC = 1.0;
1760 : // if (vTrack->GetTPCNclsF()>0) {
1761 : // ratioCrossedRowsOverFindableClustersTPC = nCrossedRowsTPC / vTrack->GetTPCNclsF();
1762 : // }
1763 :
1764 : // Int_t nClustersTPCShared = vTrack->GetTPCnclsS();
1765 : Int_t nClustersTPCShared = -1;
1766 :
1767 : // Float_t fracClustersTPCShared = -1.;
1768 :
1769 : Float_t chi2PerClusterITS = -1;
1770 : Float_t chi2PerClusterTPC = -1;
1771 : // if (nClustersITS!=0)
1772 : // chi2PerClusterITS = vTrack->GetITSchi2()/Float_t(nClustersITS);
1773 : // if (nClustersTPC!=0) {
1774 : // if(fCutRequireTPCStandAlone) {
1775 : // chi2PerClusterTPC = vTrack->GetTPCchi2Iter1()/Float_t(nClustersTPC);
1776 : // } else {
1777 : // chi2PerClusterTPC = vTrack->GetTPCchi2()/Float_t(nClustersTPC);
1778 : // }
1779 : // // fracClustersTPCShared = Float_t(nClustersTPCShared)/Float_t(nClustersTPC);
1780 : // }
1781 :
1782 0 : Double_t extCov[15];
1783 0 : for(int i=0; i<15; i++) extCov[i]=0.;
1784 : // vTrack->GetExternalCovariance(extCov);
1785 :
1786 0 : Float_t b[2];
1787 0 : Float_t bCov[3];
1788 0 : vTrack->GetImpactParameters(b,bCov);
1789 0 : if (bCov[0]<=0 || bCov[2]<=0) {
1790 0 : AliDebug(1, "Estimated b resolution lower or equal zero!");
1791 0 : bCov[0]=0; bCov[2]=0;
1792 0 : }
1793 :
1794 :
1795 : // set pt-dependent DCA cuts, if requested
1796 0 : SetPtDepDCACuts(vTrack->Pt());
1797 :
1798 :
1799 0 : Float_t dcaToVertexXY = b[0];
1800 0 : Float_t dcaToVertexZ = b[1];
1801 :
1802 : Float_t dcaToVertex = -1;
1803 :
1804 0 : if (fCutDCAToVertex2D)
1805 : {
1806 0 : dcaToVertex = TMath::Sqrt(dcaToVertexXY*dcaToVertexXY/fCutMaxDCAToVertexXY/fCutMaxDCAToVertexXY + dcaToVertexZ*dcaToVertexZ/fCutMaxDCAToVertexZ/fCutMaxDCAToVertexZ);
1807 0 : }
1808 : else
1809 0 : dcaToVertex = TMath::Sqrt(dcaToVertexXY*dcaToVertexXY + dcaToVertexZ*dcaToVertexZ);
1810 :
1811 : // getting the kinematic variables of the track
1812 : // (assuming the mass is known)
1813 0 : Double_t p[3];
1814 : // vTrack->GetPxPyPz(p);
1815 0 : vTrack->PxPyPz(p);
1816 :
1817 : // Changed from float to double to prevent rounding errors leading to negative
1818 : // log arguments (M.G.)
1819 0 : Double_t momentum = TMath::Sqrt(p[0]*p[0] + p[1]*p[1] + p[2]*p[2]);
1820 0 : Double_t pt = TMath::Sqrt(p[0]*p[0] + p[1]*p[1]);
1821 : // Double_t mass = vTrack->GetMass();
1822 : // Double_t energy = TMath::Sqrt(mass*mass + momentum*momentum);
1823 :
1824 : //y-eta related calculations
1825 : // Float_t eta = -100.;
1826 : // Float_t y = -100.;
1827 : // if((momentum != TMath::Abs(p[2]))&&(momentum != 0))
1828 : // eta = 0.5*TMath::Log((momentum + p[2])/(momentum - p[2]));
1829 : // if((energy != TMath::Abs(p[2]))&&(energy != 0))
1830 : // y = 0.5*TMath::Log((energy + p[2])/(energy - p[2]));
1831 0 : Float_t eta = vTrack->Eta();
1832 0 : Float_t y = vTrack->Y();
1833 :
1834 : // if (extCov[14] < 0)
1835 : // {
1836 : // AliWarning(Form("GetSigma1Pt2() returns negative value for external covariance matrix element fC[14]: %f. Corrupted track information, track will not be accepted!", extCov[14]));
1837 : // return kFALSE;
1838 : // }
1839 0 : Float_t relUncertainty1Pt = TMath::Sqrt(extCov[14])*pt;
1840 :
1841 : //########################################################################
1842 : // cut the track?
1843 :
1844 0 : Bool_t cuts[kNCuts];
1845 0 : for (Int_t i=0; i<kNCuts; i++) cuts[i]=kFALSE;
1846 :
1847 : // track quality cuts
1848 0 : if (fCutRequireTPCRefit && (status&AliESDtrack::kTPCrefit)==0)
1849 0 : cuts[0]=kTRUE;
1850 : // if (fCutRequireTPCStandAlone && (status&AliESDtrack::kTPCin)==0)
1851 : // cuts[1]=kTRUE;
1852 0 : if (fCutRequireITSRefit && (status&AliESDtrack::kITSrefit)==0)
1853 0 : cuts[2]=kTRUE;
1854 0 : if (nClustersTPC<fCutMinNClusterTPC)
1855 0 : cuts[3]=kTRUE;
1856 0 : if (nClustersITS<fCutMinNClusterITS)
1857 0 : cuts[4]=kTRUE;
1858 : // if (chi2PerClusterTPC>fCutMaxChi2PerClusterTPC)
1859 : // cuts[5]=kTRUE;
1860 : // if (chi2PerClusterITS>fCutMaxChi2PerClusterITS)
1861 : // cuts[6]=kTRUE;
1862 : // if (extCov[0] > fCutMaxC11)
1863 : // cuts[7]=kTRUE;
1864 : // if (extCov[2] > fCutMaxC22)
1865 : // cuts[8]=kTRUE;
1866 : // if (extCov[5] > fCutMaxC33)
1867 : // cuts[9]=kTRUE;
1868 : // if (extCov[9] > fCutMaxC44)
1869 : // cuts[10]=kTRUE;
1870 : // if (extCov[14] > fCutMaxC55)
1871 : // cuts[11]=kTRUE;
1872 :
1873 : // cut 12 and 13 see below
1874 :
1875 : // if (!fCutAcceptKinkDaughters && vTrack->GetKinkIndex(0)>0)
1876 : // cuts[14]=kTRUE;
1877 : // track kinematics cut
1878 0 : if((momentum < fPMin) || (momentum > fPMax))
1879 0 : cuts[15]=kTRUE;
1880 0 : if((pt < fPtMin) || (pt > fPtMax))
1881 0 : cuts[16] = kTRUE;
1882 0 : if((p[0] < fPxMin) || (p[0] > fPxMax))
1883 0 : cuts[17] = kTRUE;
1884 0 : if((p[1] < fPyMin) || (p[1] > fPyMax))
1885 0 : cuts[18] = kTRUE;
1886 0 : if((p[2] < fPzMin) || (p[2] > fPzMax))
1887 0 : cuts[19] = kTRUE;
1888 0 : if((eta < fEtaMin) || (eta > fEtaMax))
1889 0 : cuts[20] = kTRUE;
1890 0 : if((y < fRapMin) || (y > fRapMax))
1891 0 : cuts[21] = kTRUE;
1892 0 : if (fCutDCAToVertex2D && dcaToVertex > 1)
1893 0 : cuts[22] = kTRUE;
1894 0 : if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexXY) > fCutMaxDCAToVertexXY)
1895 0 : cuts[23] = kTRUE;
1896 0 : if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexZ) > fCutMaxDCAToVertexZ)
1897 0 : cuts[24] = kTRUE;
1898 0 : if (fCutDCAToVertex2D && fCutMinDCAToVertexXY > 0 && fCutMinDCAToVertexZ > 0 && dcaToVertexXY*dcaToVertexXY/fCutMinDCAToVertexXY/fCutMinDCAToVertexXY + dcaToVertexZ*dcaToVertexZ/fCutMinDCAToVertexZ/fCutMinDCAToVertexZ < 1)
1899 0 : cuts[25] = kTRUE;
1900 0 : if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexXY) < fCutMinDCAToVertexXY)
1901 0 : cuts[26] = kTRUE;
1902 0 : if (!fCutDCAToVertex2D && TMath::Abs(dcaToVertexZ) < fCutMinDCAToVertexZ)
1903 0 : cuts[27] = kTRUE;
1904 :
1905 : // for (Int_t i = 0; i < 3; i++) {
1906 : // if(!(vTrack->GetStatus()&AliESDtrack::kITSupg)) { // current ITS
1907 : // cuts[28+i] = !CheckITSClusterRequirement(fCutClusterRequirementITS[i], vTrack->HasPointOnITSLayer(i*2), vTrack->HasPointOnITSLayer(i*2+1));
1908 : // } else { // upgraded ITS (7 layers)
1909 : // // at the moment, for L012 the layers 12 are considered together
1910 : // if(i==0) { // L012
1911 : // cuts[28+i] = !CheckITSClusterRequirement(fCutClusterRequirementITS[i], vTrack->HasPointOnITSLayer(0), (vTrack->HasPointOnITSLayer(1))&(vTrack->HasPointOnITSLayer(2)));
1912 : // } else { // L34 or L56
1913 : // cuts[28+i] = !CheckITSClusterRequirement(fCutClusterRequirementITS[i], vTrack->HasPointOnITSLayer(i*2+1), vTrack->HasPointOnITSLayer(i*2+2));
1914 : // }
1915 : // }
1916 : // }
1917 :
1918 0 : if(fCutRequireITSStandAlone || fCutRequireITSpureSA){
1919 0 : if ((status & AliESDtrack::kITSin) == 0 || (status & AliESDtrack::kTPCin)){
1920 : // TPC tracks
1921 0 : cuts[31] = kTRUE;
1922 0 : }else{
1923 : // ITS standalone tracks
1924 0 : if(fCutRequireITSStandAlone && !fCutRequireITSpureSA){
1925 0 : if(status & AliESDtrack::kITSpureSA) cuts[31] = kTRUE;
1926 0 : }else if(fCutRequireITSpureSA){
1927 0 : if(!(status & AliESDtrack::kITSpureSA)) cuts[31] = kTRUE;
1928 : }
1929 : }
1930 : }
1931 :
1932 : // if (relUncertainty1Pt > fCutMaxRel1PtUncertainty)
1933 : // cuts[32] = kTRUE;
1934 :
1935 : // if (!fCutAcceptSharedTPCClusters && nClustersTPCShared!=0)
1936 : // cuts[33] = kTRUE;
1937 :
1938 : // if (fracClustersTPCShared > fCutMaxFractionSharedTPCClusters)
1939 : // cuts[34] = kTRUE;
1940 :
1941 : Int_t nITSPointsForPid=0;
1942 0 : UChar_t clumap=vTrack->GetITSClusterMap();
1943 0 : for(Int_t i=2; i<6; i++){
1944 0 : if(clumap&(1<<i)) ++nITSPointsForPid;
1945 : }
1946 0 : if(fCutRequireITSPid && nITSPointsForPid<3) cuts[35] = kTRUE;
1947 :
1948 :
1949 0 : if (nCrossedRowsTPC<fCutMinNCrossedRowsTPC)
1950 0 : cuts[36]=kTRUE;
1951 : // if (ratioCrossedRowsOverFindableClustersTPC<fCutMinRatioCrossedRowsOverFindableClustersTPC)
1952 : // cuts[37]=kTRUE;
1953 :
1954 : Int_t nMissITSpts=0;
1955 : // Int_t idet,statusLay;
1956 : // Float_t xloc,zloc;
1957 : // for(Int_t iLay=0; iLay<6; iLay++){
1958 : // Bool_t retc=vTrack->GetITSModuleIndexInfo(iLay,idet,statusLay,xloc,zloc);
1959 : // if(retc && statusLay==5) ++nMissITSpts;
1960 : // }
1961 : // if(nMissITSpts>fCutMaxMissingITSPoints) cuts[38] = kTRUE;
1962 :
1963 : //kTOFout
1964 : // if (fCutRequireTOFout && (status&AliESDtrack::kTOFout)==0)
1965 : // cuts[40]=kTRUE;
1966 :
1967 : // TOF signal Dz cut
1968 : // Float_t dxTOF = vTrack->GetTOFsignalDx();
1969 : // Float_t dzTOF = vTrack->GetTOFsignalDz();
1970 : Float_t dxTOF = 0.;
1971 : Float_t dzTOF = 0.;
1972 : // if (fFlagCutTOFdistance && (vTrack->GetStatus() & AliESDtrack::kTOFout) == AliESDtrack::kTOFout){ // applying the TOF distance cut only if requested, and only on tracks that reached the TOF and where associated with a TOF hit
1973 : // if (fgBeamTypeFlag < 0) { // the check on the beam type was not done yet
1974 : // const AliESDEvent* event = vTrack->GetESDEvent();
1975 : // if (event){
1976 : // TString beamTypeESD = event->GetBeamType();
1977 : // AliDebug(2,Form("Beam type from ESD event = %s",beamTypeESD.Data()));
1978 : // if (beamTypeESD.CompareTo("A-A",TString::kIgnoreCase) == 0){ // we are in PbPb collisions --> fgBeamTypeFlag will be set to 1, to apply the cut on TOF signal Dz
1979 : // fgBeamTypeFlag = 1;
1980 : // }
1981 : // else { // we are NOT in PbPb collisions --> fgBeamTypeFlag will be set to 0, to NOT apply the cu6 on TOF signal Dz
1982 : // fgBeamTypeFlag = 0;
1983 : // }
1984 : // }
1985 : // else{
1986 : // AliFatal("Beam type not available, but it is needed to apply the TOF cut!");
1987 : // }
1988 : // }
1989 :
1990 : // if (fgBeamTypeFlag == 1){ // we are in PbPb collisions --> apply the cut on TOF signal Dz
1991 : // Float_t radiusTOF = TMath::Sqrt(dxTOF*dxTOF + dzTOF*dzTOF);
1992 : // AliDebug(3,Form("TOF check (with fCutTOFdistance = %f) --> dx = %f, dz = %f, radius = %f", fCutTOFdistance, dxTOF, dzTOF, radiusTOF));
1993 : // if (radiusTOF > fCutTOFdistance){
1994 : // AliDebug(2, Form("************* the radius is outside the range! %f > %f, the track will be skipped", radiusTOF, fCutTOFdistance));
1995 : // cuts[41] = kTRUE;
1996 : // }
1997 : // }
1998 : // }
1999 :
2000 : Bool_t cut=kFALSE;
2001 0 : for (Int_t i=0; i<kNCuts; i++)
2002 0 : if (cuts[i]) {cut = kTRUE;}
2003 :
2004 : // for performance evaluate the CPU intensive cuts only when the others have passed, and when they are requested
2005 : Double_t chi2TPCConstrainedVsGlobal = -2;
2006 : Float_t nSigmaToVertex = -2;
2007 0 : if (!cut)
2008 : {
2009 : // getting the track to vertex parameters
2010 0 : if (fCutSigmaToVertexRequired)
2011 : {
2012 0 : nSigmaToVertex = GetSigmaToVertexVTrack(vTrack);
2013 0 : if (nSigmaToVertex > fCutNsigmaToVertex && fCutSigmaToVertexRequired)
2014 : {
2015 0 : cuts[12] = kTRUE;
2016 : cut = kTRUE;
2017 0 : }
2018 : // if n sigma could not be calculated
2019 0 : if (nSigmaToVertex<0 && fCutSigmaToVertexRequired)
2020 : {
2021 0 : cuts[13] = kTRUE;
2022 : cut = kTRUE;
2023 0 : }
2024 : }
2025 :
2026 : // // max chi2 TPC constrained vs global track only if track passed the other cut
2027 : // if (fCutMaxChi2TPCConstrainedVsGlobal < 1e9)
2028 : // {
2029 : // const AliESDEvent* esdEvent = vTrack->GetESDEvent();
2030 :
2031 : // if (!esdEvent)
2032 : // AliFatal("fCutMaxChi2TPCConstrainedVsGlobal set but ESD event not set in AliESDTrack. Use AliESDTrack::SetESDEvent before calling AliESDtrackCuts.");
2033 :
2034 : // // get vertex
2035 : // const AliESDVertex* vertex = 0;
2036 : // if (fCutMaxChi2TPCConstrainedVsGlobalVertexType & kVertexTracks)
2037 : // vertex = esdEvent->GetPrimaryVertexTracks();
2038 :
2039 : // if ((!vertex || !vertex->GetStatus()) && fCutMaxChi2TPCConstrainedVsGlobalVertexType & kVertexSPD)
2040 : // vertex = esdEvent->GetPrimaryVertexSPD();
2041 :
2042 : // if ((!vertex || !vertex->GetStatus()) && fCutMaxChi2TPCConstrainedVsGlobalVertexType & kVertexTPC)
2043 : // vertex = esdEvent->GetPrimaryVertexTPC();
2044 :
2045 : // if (vertex->GetStatus())
2046 : // chi2TPCConstrainedVsGlobal = vTrack->GetChi2TPCConstrainedVsGlobal(vertex);
2047 :
2048 : // if (chi2TPCConstrainedVsGlobal < 0 || chi2TPCConstrainedVsGlobal > fCutMaxChi2TPCConstrainedVsGlobal)
2049 : // {
2050 : // cuts[39] = kTRUE;
2051 : // cut = kTRUE;
2052 : // }
2053 : // }
2054 :
2055 : // // max length in active volume
2056 : // Float_t lengthInActiveZoneTPC = -1;
2057 : // if (fCutMinLengthActiveVolumeTPC > 1.) { // do the calculation only if needed to save cpu-time
2058 : // if (vTrack->GetESDEvent()) {
2059 : // if (vTrack->GetInnerParam()) lengthInActiveZoneTPC = vTrack->GetLengthInActiveZone(1, 1.8, 220, vTrack->GetESDEvent()->GetMagneticField());
2060 : // //
2061 : // if (lengthInActiveZoneTPC < fCutMinLengthActiveVolumeTPC ) {
2062 : // cuts[42] = kTRUE;
2063 : // cut = kTRUE;
2064 : // }
2065 : // }
2066 : // }
2067 :
2068 :
2069 : }
2070 :
2071 : //########################################################################
2072 : // filling histograms
2073 0 : if (fHistogramsOn) {
2074 0 : fhCutStatistics->Fill(fhCutStatistics->GetBinCenter(fhCutStatistics->GetXaxis()->FindBin("n tracks")));
2075 0 : if (cut)
2076 0 : fhCutStatistics->Fill(fhCutStatistics->GetBinCenter(fhCutStatistics->GetXaxis()->FindBin("n cut tracks")));
2077 :
2078 0 : for (Int_t i=0; i<kNCuts; i++) {
2079 0 : if (fhCutStatistics->GetXaxis()->FindBin(fgkCutNames[i]) < 1)
2080 0 : AliFatal(Form("Inconsistency! Cut %d with name %s not found", i, fgkCutNames[i]));
2081 :
2082 0 : if (cuts[i])
2083 0 : fhCutStatistics->Fill(fhCutStatistics->GetBinCenter(fhCutStatistics->GetXaxis()->FindBin(fgkCutNames[i])));
2084 :
2085 0 : for (Int_t j=i; j<kNCuts; j++) {
2086 0 : if (cuts[i] && cuts[j]) {
2087 0 : Float_t xC = fhCutCorrelation->GetXaxis()->GetBinCenter(fhCutCorrelation->GetXaxis()->FindBin(fgkCutNames[i]));
2088 0 : Float_t yC = fhCutCorrelation->GetYaxis()->GetBinCenter(fhCutCorrelation->GetYaxis()->FindBin(fgkCutNames[j]));
2089 0 : fhCutCorrelation->Fill(xC, yC);
2090 0 : }
2091 : }
2092 : }
2093 0 : }
2094 :
2095 : // now we loop over the filling of the histograms twice: once "before" the cut, once "after"
2096 : // the code is not in a function due to too many local variables that would need to be passed
2097 :
2098 0 : for (Int_t id = 0; id < 2; id++)
2099 : {
2100 : // id = 0 --> before cut
2101 : // id = 1 --> after cut
2102 :
2103 0 : if (fHistogramsOn)
2104 : {
2105 0 : fhNClustersITS[id]->Fill(nClustersITS);
2106 0 : fhNClustersTPC[id]->Fill(nClustersTPC);
2107 0 : fhNSharedClustersTPC[id]->Fill(nClustersTPCShared);
2108 0 : fhNCrossedRowsTPC[id]->Fill(nCrossedRowsTPC);
2109 0 : fhRatioCrossedRowsOverFindableClustersTPC[id]->Fill(ratioCrossedRowsOverFindableClustersTPC);
2110 0 : fhChi2PerClusterITS[id]->Fill(chi2PerClusterITS);
2111 0 : fhChi2PerClusterTPC[id]->Fill(chi2PerClusterTPC);
2112 0 : fhChi2TPCConstrainedVsGlobal[id]->Fill(chi2TPCConstrainedVsGlobal);
2113 0 : fhNClustersForITSPID[id]->Fill(nITSPointsForPid);
2114 0 : fhNMissingITSPoints[id]->Fill(nMissITSpts);
2115 :
2116 0 : fhC11[id]->Fill(extCov[0]);
2117 0 : fhC22[id]->Fill(extCov[2]);
2118 0 : fhC33[id]->Fill(extCov[5]);
2119 0 : fhC44[id]->Fill(extCov[9]);
2120 0 : fhC55[id]->Fill(extCov[14]);
2121 :
2122 0 : fhRel1PtUncertainty[id]->Fill(relUncertainty1Pt);
2123 :
2124 0 : fhPt[id]->Fill(pt);
2125 0 : fhEta[id]->Fill(eta);
2126 0 : fhTOFdistance[id]->Fill(dxTOF, dzTOF);
2127 :
2128 : Float_t bRes[2];
2129 0 : bRes[0] = TMath::Sqrt(bCov[0]);
2130 0 : bRes[1] = TMath::Sqrt(bCov[2]);
2131 :
2132 0 : fhDZ[id]->Fill(b[1]);
2133 0 : fhDXY[id]->Fill(b[0]);
2134 0 : fhDXYDZ[id]->Fill(dcaToVertex);
2135 0 : fhDXYvsDZ[id]->Fill(b[1],b[0]);
2136 :
2137 0 : if (bRes[0]!=0 && bRes[1]!=0) {
2138 0 : fhDZNormalized[id]->Fill(b[1]/bRes[1]);
2139 0 : fhDXYNormalized[id]->Fill(b[0]/bRes[0]);
2140 0 : fhDXYvsDZNormalized[id]->Fill(b[1]/bRes[1], b[0]/bRes[0]);
2141 0 : fhNSigmaToVertex[id]->Fill(nSigmaToVertex);
2142 0 : }
2143 0 : }
2144 :
2145 : // cut the track
2146 0 : if (cut)
2147 0 : return kFALSE;
2148 : }
2149 :
2150 0 : return kTRUE;
2151 0 : }
2152 :
2153 :
2154 :
2155 :
2156 :
2157 :
2158 :
2159 : //____________________________________________________________________
2160 : Bool_t AliESDtrackCuts::CheckITSClusterRequirement(ITSClusterRequirement req, Bool_t clusterL1, Bool_t clusterL2)
2161 : {
2162 : // checks if the cluster requirement is fullfilled (in this case: return kTRUE)
2163 :
2164 3132 : switch (req)
2165 : {
2166 2838 : case kOff: return kTRUE;
2167 54 : case kNone: return !clusterL1 && !clusterL2;
2168 276 : case kAny: return clusterL1 || clusterL2;
2169 0 : case kFirst: return clusterL1;
2170 0 : case kOnlyFirst: return clusterL1 && !clusterL2;
2171 0 : case kSecond: return clusterL2;
2172 0 : case kOnlySecond: return clusterL2 && !clusterL1;
2173 0 : case kBoth: return clusterL1 && clusterL2;
2174 : }
2175 :
2176 0 : return kFALSE;
2177 3132 : }
2178 :
2179 : //____________________________________________________________________
2180 : AliESDtrack* AliESDtrackCuts::GetTPCOnlyTrack(const AliESDEvent* esd, Int_t iTrack)
2181 : {
2182 : // Utility function to create a TPC only track from the given esd track
2183 : //
2184 : // IMPORTANT: The track has to be deleted by the user
2185 : //
2186 : // NB. most of the functionality to get a TPC only track from an ESD track is in AliESDtrack, where it should be
2187 : // there are only missing propagations here that are needed for old data
2188 : // this function will therefore become obsolete
2189 : //
2190 : // adapted from code provided by CKB
2191 :
2192 0 : if (!esd->GetPrimaryVertexTPC())
2193 0 : return 0; // No TPC vertex no TPC tracks
2194 :
2195 0 : if(!esd->GetPrimaryVertexTPC()->GetStatus())
2196 0 : return 0; // TPC Vertex is created by default in AliESDEvent, do not use in this case
2197 :
2198 0 : AliESDtrack* track = esd->GetTrack(iTrack);
2199 0 : if (!track)
2200 0 : return 0;
2201 :
2202 0 : AliESDtrack *tpcTrack = new AliESDtrack();
2203 :
2204 : // only true if we have a tpc track
2205 0 : if (!track->FillTPCOnlyTrack(*tpcTrack))
2206 : {
2207 0 : delete tpcTrack;
2208 0 : return 0;
2209 : }
2210 :
2211 0 : return tpcTrack;
2212 0 : }
2213 :
2214 : //____________________________________________________________________
2215 : AliESDtrack* AliESDtrackCuts::GetTPCOnlyTrackFromVEvent(const AliVEvent* vEvent, Int_t iTrack)
2216 : {
2217 : // Utility function to create a TPC only track from the given track in a VEvent. This will return 0 if the VEvent isnt an ESDEvent
2218 : //
2219 : // IMPORTANT: The track has to be deleted by the user
2220 : //
2221 : // NB. most of the functionality to get a TPC only track from an ESD track is in AliESDtrack, where it should be
2222 : // there are only missing propagations here that are needed for old data
2223 : // this function will therefore become obsolete
2224 : //
2225 : // adapted from code provided by CKB
2226 :
2227 0 : if (!vEvent->GetPrimaryVertexTPC())
2228 0 : return 0; // No TPC vertex no TPC tracks
2229 :
2230 0 : if(!vEvent->GetPrimaryVertexTPC()->GetStatus())
2231 0 : return 0; // TPC Vertex is created by default in AliESDEvent, do not use in this case
2232 :
2233 0 : AliVParticle* vParticle = vEvent->GetTrack(iTrack);
2234 0 : if (!vParticle) return 0;
2235 :
2236 0 : AliESDtrack *esdTrack = dynamic_cast<AliESDtrack*>(vParticle);
2237 0 : if(!esdTrack) return 0;
2238 :
2239 0 : AliESDtrack *tpcTrack = new AliESDtrack();
2240 :
2241 : // only true if we have a tpc track
2242 0 : if (!esdTrack->FillTPCOnlyTrack(*tpcTrack))
2243 : {
2244 0 : delete tpcTrack;
2245 0 : return 0;
2246 : }
2247 :
2248 0 : return tpcTrack;
2249 0 : }
2250 :
2251 : //____________________________________________________________________
2252 : TObjArray* AliESDtrackCuts::GetAcceptedTracks(const AliESDEvent* esd, Bool_t bTPC)
2253 : {
2254 : //
2255 : // returns an array of all tracks that pass the cuts
2256 : // or an array of TPC only tracks (propagated to the TPC vertex during reco)
2257 : // tracks that pass the cut
2258 : //
2259 : // NOTE: List has to be deleted by the user
2260 :
2261 0 : TObjArray* acceptedTracks = new TObjArray();
2262 :
2263 : // loop over esd tracks
2264 0 : for (Int_t iTrack = 0; iTrack < esd->GetNumberOfTracks(); iTrack++) {
2265 0 : if(bTPC){
2266 0 : if(!esd->GetPrimaryVertexTPC())return acceptedTracks; // No TPC vertex no TPC tracks
2267 0 : if(!esd->GetPrimaryVertexTPC()->GetStatus())return acceptedTracks; // No proper TPC vertex, only the default
2268 :
2269 0 : AliESDtrack *tpcTrack = GetTPCOnlyTrack(esd, iTrack);
2270 0 : if (!tpcTrack)
2271 0 : continue;
2272 :
2273 0 : if (AcceptTrack(tpcTrack)) {
2274 0 : acceptedTracks->Add(tpcTrack);
2275 0 : }
2276 : else
2277 0 : delete tpcTrack;
2278 0 : }
2279 : else
2280 : {
2281 0 : AliESDtrack* track = esd->GetTrack(iTrack);
2282 0 : if(AcceptTrack(track))
2283 0 : acceptedTracks->Add(track);
2284 : }
2285 : }
2286 0 : if(bTPC)acceptedTracks->SetOwner(kTRUE);
2287 0 : return acceptedTracks;
2288 0 : }
2289 :
2290 : //____________________________________________________________________
2291 : Int_t AliESDtrackCuts::CountAcceptedTracks(const AliESDEvent* const esd)
2292 : {
2293 : //
2294 : // returns an the number of tracks that pass the cuts
2295 : //
2296 :
2297 : Int_t count = 0;
2298 :
2299 : // loop over esd tracks
2300 308 : for (Int_t iTrack = 0; iTrack < esd->GetNumberOfTracks(); iTrack++) {
2301 142 : AliESDtrack* track = esd->GetTrack(iTrack);
2302 142 : if (AcceptTrack(track))
2303 96 : count++;
2304 : }
2305 :
2306 8 : return count;
2307 : }
2308 :
2309 : //____________________________________________________________________
2310 : void AliESDtrackCuts::DefineHistograms(Int_t color) {
2311 : //
2312 : // diagnostics histograms are defined
2313 : //
2314 :
2315 0 : fHistogramsOn=kTRUE;
2316 :
2317 0 : Bool_t oldStatus = TH1::AddDirectoryStatus();
2318 0 : TH1::AddDirectory(kFALSE);
2319 :
2320 : //###################################################################################
2321 : // defining histograms
2322 :
2323 0 : fhCutStatistics = new TH1F("cut_statistics","cut statistics",kNCuts+4,-0.5,kNCuts+3.5);
2324 :
2325 0 : fhCutStatistics->GetXaxis()->SetBinLabel(1,"n tracks");
2326 0 : fhCutStatistics->GetXaxis()->SetBinLabel(2,"n cut tracks");
2327 :
2328 0 : fhCutCorrelation = new TH2F("cut_correlation","cut correlation",kNCuts,-0.5,kNCuts-0.5,kNCuts,-0.5,kNCuts-0.5);;
2329 :
2330 0 : for (Int_t i=0; i<kNCuts; i++) {
2331 0 : fhCutStatistics->GetXaxis()->SetBinLabel(i+4,fgkCutNames[i]);
2332 0 : fhCutCorrelation->GetXaxis()->SetBinLabel(i+1,fgkCutNames[i]);
2333 0 : fhCutCorrelation->GetYaxis()->SetBinLabel(i+1,fgkCutNames[i]);
2334 : }
2335 :
2336 0 : fhCutStatistics ->SetLineColor(color);
2337 0 : fhCutCorrelation ->SetLineColor(color);
2338 0 : fhCutStatistics ->SetLineWidth(2);
2339 0 : fhCutCorrelation ->SetLineWidth(2);
2340 :
2341 0 : for (Int_t i=0; i<2; i++) {
2342 0 : fhNClustersITS[i] = new TH1F("nClustersITS" ,"",8,-0.5,7.5);
2343 0 : fhNClustersTPC[i] = new TH1F("nClustersTPC" ,"",165,-0.5,164.5);
2344 0 : fhNSharedClustersTPC[i] = new TH1F("nSharedClustersTPC" ,"",165,-0.5,164.5);
2345 0 : fhNCrossedRowsTPC[i] = new TH1F("nCrossedRowsTPC" ,"",165,-0.5,164.5);
2346 0 : fhRatioCrossedRowsOverFindableClustersTPC[i] = new TH1F("ratioCrossedRowsOverFindableClustersTPC" ,"",60,0,1.5);
2347 0 : fhChi2PerClusterITS[i] = new TH1F("chi2PerClusterITS","",500,0,10);
2348 0 : fhChi2PerClusterTPC[i] = new TH1F("chi2PerClusterTPC","",500,0,10);
2349 0 : fhChi2TPCConstrainedVsGlobal[i] = new TH1F("chi2TPCConstrainedVsGlobal","",600,-2,50);
2350 0 : fhNClustersForITSPID[i] = new TH1F("nPointsForITSpid","",5,-0.5,4.5);
2351 0 : fhNMissingITSPoints[i] = new TH1F("nMissingITSClusters","",7,-0.5,6.5);
2352 :
2353 0 : fhC11[i] = new TH1F("covMatrixDiagonal11","",2000,0,20);
2354 0 : fhC22[i] = new TH1F("covMatrixDiagonal22","",2000,0,20);
2355 0 : fhC33[i] = new TH1F("covMatrixDiagonal33","",1000,0,0.1);
2356 0 : fhC44[i] = new TH1F("covMatrixDiagonal44","",1000,0,0.1);
2357 0 : fhC55[i] = new TH1F("covMatrixDiagonal55","",1000,0,5);
2358 :
2359 0 : fhRel1PtUncertainty[i] = new TH1F("rel1PtUncertainty","",1000,0,5);
2360 :
2361 0 : fhDXY[i] = new TH1F("dXY" ,"",500,-10,10);
2362 0 : fhDZ[i] = new TH1F("dZ" ,"",500,-10,10);
2363 0 : fhDXYDZ[i] = new TH1F("dXYDZ" ,"",500,0,10);
2364 0 : fhDXYvsDZ[i] = new TH2F("dXYvsDZ","",200,-10,10,200,-10,10);
2365 :
2366 0 : fhDXYNormalized[i] = new TH1F("dXYNormalized" ,"",500,-10,10);
2367 0 : fhDZNormalized[i] = new TH1F("dZNormalized" ,"",500,-10,10);
2368 0 : fhDXYvsDZNormalized[i] = new TH2F("dXYvsDZNormalized","",200,-10,10,200,-10,10);
2369 :
2370 0 : fhNSigmaToVertex[i] = new TH1F("nSigmaToVertex","",500,0,10);
2371 :
2372 0 : fhPt[i] = new TH1F("pt" ,"p_{T} distribution;p_{T} (GeV/c)", 800, 0.0, 10.0);
2373 0 : fhEta[i] = new TH1F("eta" ,"#eta distribution;#eta",40,-2.0,2.0);
2374 0 : fhTOFdistance[i] = new TH2F("TOFdistance" ,"TOF distance;dx (cm};dz (cm)", 150, -15, 15, 150, -15, 15);
2375 :
2376 0 : fhNClustersITS[i]->SetTitle("n ITS clusters");
2377 0 : fhNClustersTPC[i]->SetTitle("n TPC clusters");
2378 0 : fhNSharedClustersTPC[i]->SetTitle("n TPC shared clusters");
2379 0 : fhChi2PerClusterITS[i]->SetTitle("#Chi^{2} per ITS cluster");
2380 0 : fhChi2PerClusterTPC[i]->SetTitle("#Chi^{2} per TPC cluster");
2381 0 : fhChi2TPCConstrainedVsGlobal[i]->SetTitle("#Chi^{2} TPC constrained track vs global track");
2382 0 : fhNClustersForITSPID[i]->SetTitle("n ITS points for PID");
2383 0 : fhNMissingITSPoints[i]->SetTitle("n ITS layers with missing cluster");
2384 :
2385 0 : fhC11[i]->SetTitle("cov 11 : #sigma_{y}^{2} [cm^{2}]");
2386 0 : fhC22[i]->SetTitle("cov 22 : #sigma_{z}^{2} [cm^{2}]");
2387 0 : fhC33[i]->SetTitle("cov 33 : #sigma_{sin(#phi)}^{2}");
2388 0 : fhC44[i]->SetTitle("cov 44 : #sigma_{tan(#theta_{dip})}^{2}");
2389 0 : fhC55[i]->SetTitle("cov 55 : #sigma_{1/p_{T}}^{2} [(c/GeV)^2]");
2390 :
2391 0 : fhRel1PtUncertainty[i]->SetTitle("rel. uncertainty of 1/p_{T}");
2392 :
2393 0 : fhDXY[i]->SetXTitle("transverse impact parameter (cm)");
2394 0 : fhDZ[i]->SetXTitle("longitudinal impact parameter (cm)");
2395 0 : fhDXYDZ[i]->SetTitle("absolute impact parameter;sqrt(dXY**2 + dZ**2) (cm)");
2396 0 : fhDXYvsDZ[i]->SetXTitle("longitudinal impact parameter (cm)");
2397 0 : fhDXYvsDZ[i]->SetYTitle("transverse impact parameter (cm)");
2398 :
2399 0 : fhDXYNormalized[i]->SetTitle("normalized trans impact par (n#sigma)");
2400 0 : fhDZNormalized[i]->SetTitle("normalized long impact par (n#sigma)");
2401 0 : fhDXYvsDZNormalized[i]->SetTitle("normalized long impact par (n#sigma)");
2402 0 : fhDXYvsDZNormalized[i]->SetYTitle("normalized trans impact par (n#sigma)");
2403 0 : fhNSigmaToVertex[i]->SetTitle("n #sigma to vertex");
2404 :
2405 0 : fhNClustersITS[i]->SetLineColor(color); fhNClustersITS[i]->SetLineWidth(2);
2406 0 : fhNClustersTPC[i]->SetLineColor(color); fhNClustersTPC[i]->SetLineWidth(2);
2407 0 : fhNSharedClustersTPC[i]->SetLineColor(color); fhNSharedClustersTPC[i]->SetLineWidth(2);
2408 0 : fhChi2PerClusterITS[i]->SetLineColor(color); fhChi2PerClusterITS[i]->SetLineWidth(2);
2409 0 : fhChi2PerClusterTPC[i]->SetLineColor(color); fhChi2PerClusterTPC[i]->SetLineWidth(2);
2410 0 : fhChi2TPCConstrainedVsGlobal[i]->SetLineColor(color); fhChi2TPCConstrainedVsGlobal[i]->SetLineWidth(2);
2411 0 : fhNClustersForITSPID[i]->SetLineColor(color); fhNClustersForITSPID[i]->SetLineWidth(2);
2412 0 : fhNMissingITSPoints[i]->SetLineColor(color); fhNMissingITSPoints[i]->SetLineWidth(2);
2413 :
2414 0 : fhC11[i]->SetLineColor(color); fhC11[i]->SetLineWidth(2);
2415 0 : fhC22[i]->SetLineColor(color); fhC22[i]->SetLineWidth(2);
2416 0 : fhC33[i]->SetLineColor(color); fhC33[i]->SetLineWidth(2);
2417 0 : fhC44[i]->SetLineColor(color); fhC44[i]->SetLineWidth(2);
2418 0 : fhC55[i]->SetLineColor(color); fhC55[i]->SetLineWidth(2);
2419 :
2420 0 : fhRel1PtUncertainty[i]->SetLineColor(color); fhRel1PtUncertainty[i]->SetLineWidth(2);
2421 :
2422 0 : fhDXY[i]->SetLineColor(color); fhDXY[i]->SetLineWidth(2);
2423 0 : fhDZ[i]->SetLineColor(color); fhDZ[i]->SetLineWidth(2);
2424 0 : fhDXYDZ[i]->SetLineColor(color); fhDXYDZ[i]->SetLineWidth(2);
2425 :
2426 0 : fhDXYNormalized[i]->SetLineColor(color); fhDXYNormalized[i]->SetLineWidth(2);
2427 0 : fhDZNormalized[i]->SetLineColor(color); fhDZNormalized[i]->SetLineWidth(2);
2428 0 : fhNSigmaToVertex[i]->SetLineColor(color); fhNSigmaToVertex[i]->SetLineWidth(2);
2429 : }
2430 :
2431 : // The number of sigmas to the vertex is per definition gaussian
2432 0 : ffDTheoretical = new TF1("nSigmaToVertexTheoretical","([0]/2.506628274)*exp(-(x**2)/2)",0,50);
2433 0 : ffDTheoretical->SetParameter(0,1);
2434 :
2435 0 : TH1::AddDirectory(oldStatus);
2436 0 : }
2437 :
2438 : //____________________________________________________________________
2439 : Bool_t AliESDtrackCuts::LoadHistograms(const Char_t* dir)
2440 : {
2441 : //
2442 : // loads the histograms from a file
2443 : // if dir is empty a directory with the name of this object is taken (like in SaveHistogram)
2444 : //
2445 :
2446 0 : if (!dir)
2447 0 : dir = GetName();
2448 :
2449 0 : if (!gDirectory->cd(dir))
2450 0 : return kFALSE;
2451 :
2452 0 : ffDTheoretical = dynamic_cast<TF1*> (gDirectory->Get("nSigmaToVertexTheory"));
2453 :
2454 0 : fhCutStatistics = dynamic_cast<TH1F*> (gDirectory->Get("cut_statistics"));
2455 0 : fhCutCorrelation = dynamic_cast<TH2F*> (gDirectory->Get("cut_correlation"));
2456 :
2457 0 : for (Int_t i=0; i<2; i++) {
2458 0 : if (i==0)
2459 : {
2460 0 : gDirectory->cd("before_cuts");
2461 0 : }
2462 : else
2463 0 : gDirectory->cd("after_cuts");
2464 :
2465 0 : fhNClustersITS[i] = dynamic_cast<TH1F*> (gDirectory->Get("nClustersITS" ));
2466 0 : fhNClustersTPC[i] = dynamic_cast<TH1F*> (gDirectory->Get("nClustersTPC" ));
2467 0 : fhNSharedClustersTPC[i] = dynamic_cast<TH1F*> (gDirectory->Get("nSharedClustersTPC" ));
2468 0 : fhNCrossedRowsTPC[i] = dynamic_cast<TH1F*> (gDirectory->Get("nCrossedRowsTPC" ));
2469 0 : fhRatioCrossedRowsOverFindableClustersTPC[i] = dynamic_cast<TH1F*> (gDirectory->Get("ratioCrossedRowsOverFindableClustersTPC" ));
2470 0 : fhChi2PerClusterITS[i] = dynamic_cast<TH1F*> (gDirectory->Get("chi2PerClusterITS"));
2471 0 : fhChi2PerClusterTPC[i] = dynamic_cast<TH1F*> (gDirectory->Get("chi2PerClusterTPC"));
2472 0 : fhChi2TPCConstrainedVsGlobal[i] = dynamic_cast<TH1F*> (gDirectory->Get("fhChi2TPCConstrainedVsGlobal"));
2473 0 : fhNClustersForITSPID[i] = dynamic_cast<TH1F*> (gDirectory->Get("nPointsForITSpid"));
2474 0 : fhNMissingITSPoints[i] = dynamic_cast<TH1F*> (gDirectory->Get("nMissingITSClusters"));
2475 :
2476 0 : fhC11[i] = dynamic_cast<TH1F*> (gDirectory->Get("covMatrixDiagonal11"));
2477 0 : fhC22[i] = dynamic_cast<TH1F*> (gDirectory->Get("covMatrixDiagonal22"));
2478 0 : fhC33[i] = dynamic_cast<TH1F*> (gDirectory->Get("covMatrixDiagonal33"));
2479 0 : fhC44[i] = dynamic_cast<TH1F*> (gDirectory->Get("covMatrixDiagonal44"));
2480 0 : fhC55[i] = dynamic_cast<TH1F*> (gDirectory->Get("covMatrixDiagonal55"));
2481 :
2482 0 : fhRel1PtUncertainty[i] = dynamic_cast<TH1F*> (gDirectory->Get("rel1PtUncertainty"));
2483 :
2484 0 : fhDXY[i] = dynamic_cast<TH1F*> (gDirectory->Get("dXY" ));
2485 0 : fhDZ[i] = dynamic_cast<TH1F*> (gDirectory->Get("dZ" ));
2486 0 : fhDXYDZ[i] = dynamic_cast<TH1F*> (gDirectory->Get("dXYDZ"));
2487 0 : fhDXYvsDZ[i] = dynamic_cast<TH2F*> (gDirectory->Get("dXYvsDZ"));
2488 :
2489 0 : fhDXYNormalized[i] = dynamic_cast<TH1F*> (gDirectory->Get("dXYNormalized" ));
2490 0 : fhDZNormalized[i] = dynamic_cast<TH1F*> (gDirectory->Get("dZNormalized" ));
2491 0 : fhDXYvsDZNormalized[i] = dynamic_cast<TH2F*> (gDirectory->Get("dXYvsDZNormalized"));
2492 0 : fhNSigmaToVertex[i] = dynamic_cast<TH1F*> (gDirectory->Get("nSigmaToVertex"));
2493 :
2494 0 : fhPt[i] = dynamic_cast<TH1F*> (gDirectory->Get("pt"));
2495 0 : fhEta[i] = dynamic_cast<TH1F*> (gDirectory->Get("eta"));
2496 0 : fhTOFdistance[i] = dynamic_cast<TH2F*> (gDirectory->Get("TOFdistance"));
2497 :
2498 0 : gDirectory->cd("../");
2499 : }
2500 :
2501 0 : gDirectory->cd("..");
2502 :
2503 0 : return kTRUE;
2504 0 : }
2505 :
2506 : //____________________________________________________________________
2507 : void AliESDtrackCuts::SaveHistograms(const Char_t* dir) {
2508 : //
2509 : // saves the histograms in a directory (dir)
2510 : //
2511 :
2512 0 : if (!fHistogramsOn) {
2513 0 : AliDebug(0, "Histograms not on - cannot save histograms!!!");
2514 : return;
2515 : }
2516 :
2517 0 : if (!dir)
2518 0 : dir = GetName();
2519 :
2520 0 : gDirectory->mkdir(dir);
2521 0 : gDirectory->cd(dir);
2522 :
2523 0 : gDirectory->mkdir("before_cuts");
2524 0 : gDirectory->mkdir("after_cuts");
2525 :
2526 : // a factor of 2 is needed since n sigma is positive
2527 0 : ffDTheoretical->SetParameter(0,2*fhNSigmaToVertex[0]->Integral("width"));
2528 0 : ffDTheoretical->Write("nSigmaToVertexTheory");
2529 :
2530 0 : fhCutStatistics->Write();
2531 0 : fhCutCorrelation->Write();
2532 :
2533 0 : for (Int_t i=0; i<2; i++) {
2534 0 : if (i==0)
2535 0 : gDirectory->cd("before_cuts");
2536 : else
2537 0 : gDirectory->cd("after_cuts");
2538 :
2539 0 : fhNClustersITS[i] ->Write();
2540 0 : fhNClustersTPC[i] ->Write();
2541 0 : fhNSharedClustersTPC[i] ->Write();
2542 0 : fhNCrossedRowsTPC[i] ->Write();
2543 0 : fhRatioCrossedRowsOverFindableClustersTPC[i] ->Write();
2544 0 : fhChi2PerClusterITS[i] ->Write();
2545 0 : fhChi2PerClusterTPC[i] ->Write();
2546 0 : fhChi2TPCConstrainedVsGlobal[i] ->Write();
2547 0 : fhNClustersForITSPID[i] ->Write();
2548 0 : fhNMissingITSPoints[i] ->Write();
2549 :
2550 0 : fhC11[i] ->Write();
2551 0 : fhC22[i] ->Write();
2552 0 : fhC33[i] ->Write();
2553 0 : fhC44[i] ->Write();
2554 0 : fhC55[i] ->Write();
2555 :
2556 0 : fhRel1PtUncertainty[i] ->Write();
2557 :
2558 0 : fhDXY[i] ->Write();
2559 0 : fhDZ[i] ->Write();
2560 0 : fhDXYDZ[i] ->Write();
2561 0 : fhDXYvsDZ[i] ->Write();
2562 :
2563 0 : fhDXYNormalized[i] ->Write();
2564 0 : fhDZNormalized[i] ->Write();
2565 0 : fhDXYvsDZNormalized[i] ->Write();
2566 0 : fhNSigmaToVertex[i] ->Write();
2567 :
2568 0 : fhPt[i] ->Write();
2569 0 : fhEta[i] ->Write();
2570 0 : fhTOFdistance[i] ->Write();
2571 :
2572 0 : gDirectory->cd("../");
2573 : }
2574 :
2575 0 : gDirectory->cd("../");
2576 0 : }
2577 :
2578 : //____________________________________________________________________
2579 : void AliESDtrackCuts::DrawHistograms()
2580 : {
2581 : // draws some histograms
2582 :
2583 0 : TCanvas* canvas1 = new TCanvas(Form("%s_1", GetName()), "Track Quality Results1", 800, 800);
2584 0 : canvas1->Divide(2, 2);
2585 :
2586 0 : canvas1->cd(1);
2587 0 : fhNClustersTPC[0]->SetStats(kFALSE);
2588 0 : fhNClustersTPC[0]->Draw();
2589 :
2590 0 : canvas1->cd(2);
2591 0 : fhChi2PerClusterTPC[0]->SetStats(kFALSE);
2592 0 : fhChi2PerClusterTPC[0]->Draw();
2593 :
2594 0 : canvas1->cd(3);
2595 0 : fhNSigmaToVertex[0]->SetStats(kFALSE);
2596 0 : fhNSigmaToVertex[0]->GetXaxis()->SetRangeUser(0, 10);
2597 0 : fhNSigmaToVertex[0]->Draw();
2598 :
2599 0 : canvas1->SaveAs(Form("%s_%s.gif", GetName(), canvas1->GetName()));
2600 :
2601 0 : TCanvas* canvas2 = new TCanvas(Form("%s_2", GetName()), "Track Quality Results2", 1200, 800);
2602 0 : canvas2->Divide(3, 2);
2603 :
2604 0 : canvas2->cd(1);
2605 0 : fhC11[0]->SetStats(kFALSE);
2606 0 : gPad->SetLogy();
2607 0 : fhC11[0]->Draw();
2608 :
2609 0 : canvas2->cd(2);
2610 0 : fhC22[0]->SetStats(kFALSE);
2611 0 : gPad->SetLogy();
2612 0 : fhC22[0]->Draw();
2613 :
2614 0 : canvas2->cd(3);
2615 0 : fhC33[0]->SetStats(kFALSE);
2616 0 : gPad->SetLogy();
2617 0 : fhC33[0]->Draw();
2618 :
2619 0 : canvas2->cd(4);
2620 0 : fhC44[0]->SetStats(kFALSE);
2621 0 : gPad->SetLogy();
2622 0 : fhC44[0]->Draw();
2623 :
2624 0 : canvas2->cd(5);
2625 0 : fhC55[0]->SetStats(kFALSE);
2626 0 : gPad->SetLogy();
2627 0 : fhC55[0]->Draw();
2628 :
2629 0 : canvas2->cd(6);
2630 0 : fhRel1PtUncertainty[0]->SetStats(kFALSE);
2631 0 : gPad->SetLogy();
2632 0 : fhRel1PtUncertainty[0]->Draw();
2633 :
2634 0 : canvas2->SaveAs(Form("%s_%s.gif", GetName(), canvas2->GetName()));
2635 :
2636 0 : TCanvas* canvas3 = new TCanvas(Form("%s_3", GetName()), "Track Quality Results3", 1200, 800);
2637 0 : canvas3->Divide(3, 2);
2638 :
2639 0 : canvas3->cd(1);
2640 0 : fhDXY[0]->SetStats(kFALSE);
2641 0 : gPad->SetLogy();
2642 0 : fhDXY[0]->Draw();
2643 :
2644 0 : canvas3->cd(2);
2645 0 : fhDZ[0]->SetStats(kFALSE);
2646 0 : gPad->SetLogy();
2647 0 : fhDZ[0]->Draw();
2648 :
2649 0 : canvas3->cd(3);
2650 0 : fhDXYvsDZ[0]->SetStats(kFALSE);
2651 0 : gPad->SetLogz();
2652 0 : gPad->SetRightMargin(0.15);
2653 0 : fhDXYvsDZ[0]->Draw("COLZ");
2654 :
2655 0 : canvas3->cd(4);
2656 0 : fhDXYNormalized[0]->SetStats(kFALSE);
2657 0 : gPad->SetLogy();
2658 0 : fhDXYNormalized[0]->Draw();
2659 :
2660 0 : canvas3->cd(5);
2661 0 : fhDZNormalized[0]->SetStats(kFALSE);
2662 0 : gPad->SetLogy();
2663 0 : fhDZNormalized[0]->Draw();
2664 :
2665 0 : canvas3->cd(6);
2666 0 : fhDXYvsDZNormalized[0]->SetStats(kFALSE);
2667 0 : gPad->SetLogz();
2668 0 : gPad->SetRightMargin(0.15);
2669 0 : fhDXYvsDZNormalized[0]->Draw("COLZ");
2670 :
2671 0 : canvas3->SaveAs(Form("%s_%s.gif", GetName(), canvas3->GetName()));
2672 :
2673 0 : TCanvas* canvas4 = new TCanvas(Form("%s_4", GetName()), "Track Quality Results4", 800, 500);
2674 0 : canvas4->Divide(2, 1);
2675 :
2676 0 : canvas4->cd(1);
2677 0 : fhCutStatistics->SetStats(kFALSE);
2678 0 : fhCutStatistics->LabelsOption("v");
2679 0 : gPad->SetBottomMargin(0.3);
2680 0 : fhCutStatistics->Draw();
2681 :
2682 0 : canvas4->cd(2);
2683 0 : fhCutCorrelation->SetStats(kFALSE);
2684 0 : fhCutCorrelation->LabelsOption("v");
2685 0 : gPad->SetBottomMargin(0.3);
2686 0 : gPad->SetLeftMargin(0.3);
2687 0 : fhCutCorrelation->Draw("COLZ");
2688 :
2689 0 : canvas4->SaveAs(Form("%s_%s.gif", GetName(), canvas4->GetName()));
2690 :
2691 : /*canvas->cd(1);
2692 : fhDXYvsDZNormalized[0]->SetStats(kFALSE);
2693 : fhDXYvsDZNormalized[0]->DrawCopy("COLZ");
2694 :
2695 : canvas->cd(2);
2696 : fhNClustersTPC[0]->SetStats(kFALSE);
2697 : fhNClustersTPC[0]->DrawCopy();
2698 :
2699 : canvas->cd(3);
2700 : fhChi2PerClusterITS[0]->SetStats(kFALSE);
2701 : fhChi2PerClusterITS[0]->DrawCopy();
2702 : fhChi2PerClusterITS[1]->SetLineColor(2);
2703 : fhChi2PerClusterITS[1]->DrawCopy("SAME");
2704 :
2705 : canvas->cd(4);
2706 : fhChi2PerClusterTPC[0]->SetStats(kFALSE);
2707 : fhChi2PerClusterTPC[0]->DrawCopy();
2708 : fhChi2PerClusterTPC[1]->SetLineColor(2);
2709 : fhChi2PerClusterTPC[1]->DrawCopy("SAME");*/
2710 0 : }
2711 : //--------------------------------------------------------------------------
2712 : void AliESDtrackCuts::SetPtDepDCACuts(Double_t pt) {
2713 : //
2714 : // set the pt-dependent DCA cuts
2715 : //
2716 :
2717 2088 : if(f1CutMaxDCAToVertexXYPtDep) {
2718 294 : fCutMaxDCAToVertexXY=f1CutMaxDCAToVertexXYPtDep->Eval(pt);
2719 294 : }
2720 :
2721 1044 : if(f1CutMaxDCAToVertexZPtDep) {
2722 0 : fCutMaxDCAToVertexZ=f1CutMaxDCAToVertexZPtDep->Eval(pt);
2723 0 : }
2724 :
2725 1044 : if(f1CutMinDCAToVertexXYPtDep) {
2726 0 : fCutMinDCAToVertexXY=f1CutMinDCAToVertexXYPtDep->Eval(pt);
2727 0 : }
2728 :
2729 1044 : if(f1CutMinDCAToVertexZPtDep) {
2730 0 : fCutMinDCAToVertexZ=f1CutMinDCAToVertexZPtDep->Eval(pt);
2731 0 : }
2732 :
2733 :
2734 1044 : return;
2735 : }
2736 :
2737 :
2738 :
2739 : //--------------------------------------------------------------------------
2740 : Bool_t AliESDtrackCuts::CheckPtDepDCA(TString dist,Bool_t print) const {
2741 : //
2742 : // Check the correctness of the string syntax
2743 : //
2744 : Bool_t retval=kTRUE;
2745 :
2746 4 : if(!dist.Contains("pt")) {
2747 0 : if(print) AliError("string must contain \"pt\"");
2748 : retval= kFALSE;
2749 0 : }
2750 4 : return retval;
2751 : }
2752 :
2753 : void AliESDtrackCuts::SetMaxDCAToVertexXYPtDep(const char *dist){
2754 :
2755 8 : if(f1CutMaxDCAToVertexXYPtDep){
2756 0 : delete f1CutMaxDCAToVertexXYPtDep;
2757 : // resetiing both
2758 0 : f1CutMaxDCAToVertexXYPtDep = 0;
2759 0 : fCutMaxDCAToVertexXYPtDep = "";
2760 0 : }
2761 8 : if(!CheckPtDepDCA(dist,kTRUE)){
2762 : return;
2763 : }
2764 4 : fCutMaxDCAToVertexXYPtDep = dist;
2765 4 : TString tmp(dist);
2766 4 : tmp.ReplaceAll("pt","x");
2767 16 : f1CutMaxDCAToVertexXYPtDep = new TFormula("f1CutMaxDCAToVertexXYPtDep",tmp.Data());
2768 :
2769 8 : }
2770 :
2771 : void AliESDtrackCuts::SetMaxDCAToVertexZPtDep(const char *dist){
2772 :
2773 :
2774 0 : if(f1CutMaxDCAToVertexZPtDep){
2775 0 : delete f1CutMaxDCAToVertexZPtDep;
2776 : // resetiing both
2777 0 : f1CutMaxDCAToVertexZPtDep = 0;
2778 0 : fCutMaxDCAToVertexZPtDep = "";
2779 0 : }
2780 0 : if(!CheckPtDepDCA(dist,kTRUE))return;
2781 :
2782 0 : fCutMaxDCAToVertexZPtDep = dist;
2783 0 : TString tmp(dist);
2784 0 : tmp.ReplaceAll("pt","x");
2785 0 : f1CutMaxDCAToVertexZPtDep = new TFormula("f1CutMaxDCAToVertexZPtDep",tmp.Data());
2786 :
2787 :
2788 0 : }
2789 :
2790 :
2791 : void AliESDtrackCuts::SetMinDCAToVertexXYPtDep(const char *dist){
2792 :
2793 :
2794 0 : if(f1CutMinDCAToVertexXYPtDep){
2795 0 : delete f1CutMinDCAToVertexXYPtDep;
2796 : // resetiing both
2797 0 : f1CutMinDCAToVertexXYPtDep = 0;
2798 0 : fCutMinDCAToVertexXYPtDep = "";
2799 0 : }
2800 0 : if(!CheckPtDepDCA(dist,kTRUE))return;
2801 :
2802 0 : fCutMinDCAToVertexXYPtDep = dist;
2803 0 : TString tmp(dist);
2804 0 : tmp.ReplaceAll("pt","x");
2805 0 : f1CutMinDCAToVertexXYPtDep = new TFormula("f1CutMinDCAToVertexXYPtDep",tmp.Data());
2806 :
2807 0 : }
2808 :
2809 :
2810 : void AliESDtrackCuts::SetMinDCAToVertexZPtDep(const char *dist){
2811 :
2812 :
2813 :
2814 0 : if(f1CutMinDCAToVertexZPtDep){
2815 0 : delete f1CutMinDCAToVertexZPtDep;
2816 : // resetiing both
2817 0 : f1CutMinDCAToVertexZPtDep = 0;
2818 0 : fCutMinDCAToVertexZPtDep = "";
2819 0 : }
2820 0 : if(!CheckPtDepDCA(dist,kTRUE))return;
2821 0 : fCutMinDCAToVertexZPtDep = dist;
2822 0 : TString tmp(dist);
2823 0 : tmp.ReplaceAll("pt","x");
2824 0 : f1CutMinDCAToVertexZPtDep = new TFormula("f1CutMinDCAToVertexZPtDep",tmp.Data());
2825 0 : }
2826 :
2827 : AliESDtrackCuts* AliESDtrackCuts::GetMultEstTrackCuts(MultEstTrackCuts cut)
2828 : {
2829 : // returns the multiplicity estimator track cuts objects to allow for user configuration
2830 : // upon first call the objects are created
2831 : //
2832 : // the cut defined here correspond to GetStandardITSTPCTrackCuts2010 (apart from the one for without SPD)
2833 :
2834 144 : if (!fgMultEstTrackCuts[kMultEstTrackCutGlobal])
2835 : {
2836 : // quality cut on ITS+TPC tracks
2837 4 : fgMultEstTrackCuts[kMultEstTrackCutGlobal] = new AliESDtrackCuts();
2838 2 : fgMultEstTrackCuts[kMultEstTrackCutGlobal]->SetMinNClustersTPC(70);
2839 2 : fgMultEstTrackCuts[kMultEstTrackCutGlobal]->SetMaxChi2PerClusterTPC(4);
2840 2 : fgMultEstTrackCuts[kMultEstTrackCutGlobal]->SetAcceptKinkDaughters(kFALSE);
2841 2 : fgMultEstTrackCuts[kMultEstTrackCutGlobal]->SetRequireTPCRefit(kTRUE);
2842 2 : fgMultEstTrackCuts[kMultEstTrackCutGlobal]->SetRequireITSRefit(kTRUE);
2843 : //multiplicity underestimate if we use global tracks with |eta| > 0.9
2844 2 : fgMultEstTrackCuts[kMultEstTrackCutGlobal]->SetEtaRange(-0.9, 0.9);
2845 :
2846 : // quality cut on ITS_SA tracks (complementary to ITS+TPC)
2847 4 : fgMultEstTrackCuts[kMultEstTrackCutITSSA] = new AliESDtrackCuts();
2848 2 : fgMultEstTrackCuts[kMultEstTrackCutITSSA]->SetRequireITSRefit(kTRUE);
2849 :
2850 : // primary selection for tracks with SPD hits
2851 4 : fgMultEstTrackCuts[kMultEstTrackCutDCAwSPD] = new AliESDtrackCuts();
2852 2 : fgMultEstTrackCuts[kMultEstTrackCutDCAwSPD]->SetClusterRequirementITS(AliESDtrackCuts::kSPD, AliESDtrackCuts::kAny);
2853 2 : fgMultEstTrackCuts[kMultEstTrackCutDCAwSPD]->SetMaxDCAToVertexXYPtDep("0.0182+0.0350/pt^1.01");
2854 2 : fgMultEstTrackCuts[kMultEstTrackCutDCAwSPD]->SetMaxDCAToVertexZ(2);
2855 :
2856 : // primary selection for tracks w/o SPD hits
2857 4 : fgMultEstTrackCuts[kMultEstTrackCutDCAwoSPD] = new AliESDtrackCuts();
2858 2 : fgMultEstTrackCuts[kMultEstTrackCutDCAwoSPD]->SetClusterRequirementITS(AliESDtrackCuts::kSPD, AliESDtrackCuts::kNone);
2859 2 : fgMultEstTrackCuts[kMultEstTrackCutDCAwoSPD]->SetMaxDCAToVertexXYPtDep("1.5*(0.0182+0.0350/pt^1.01)");
2860 2 : fgMultEstTrackCuts[kMultEstTrackCutDCAwoSPD]->SetMaxDCAToVertexZ(2);
2861 2 : }
2862 :
2863 72 : return fgMultEstTrackCuts[cut];
2864 0 : }
2865 :
2866 : Int_t AliESDtrackCuts::GetReferenceMultiplicity(const AliESDEvent* esd, MultEstTrackType trackType, Float_t etaRange, Float_t etaCent)
2867 : {
2868 : // Get multiplicity estimate based on TPC/ITS tracks and tracklets
2869 : // Adapted for AliESDtrackCuts from a version developed by: Ruben Shahoyan, Anton Alkin, Arvinder Palaha
2870 : //
2871 : // Returns a negative value if no reliable estimate can be provided:
2872 : // -1 SPD vertex missing
2873 : // -2 SPD VertexerZ dispersion too large
2874 : // -3 Track vertex missing (not checked for kTracklets)
2875 : // -4 Distance between SPD and track vertices too large (not checked for kTracklets)
2876 : //
2877 : // WARNING This functions does not cut on the z vtx. Depending on the eta range requested, you need to restrict your z vertex range!
2878 : //
2879 : // Strategy for combined estimators
2880 : // 1. Count global tracks and flag them
2881 : // 2. Count ITSSA as complementaries for ITSTPC+ or as main tracks
2882 : // 3. Count the complementary SPD tracklets
2883 : //
2884 : // Estimator is calculated for etaCent-etaRange : etaCent+etaRange
2885 : //
2886 : const AliESDVertex* vertices[2];
2887 48 : vertices[0] = esd->GetPrimaryVertexSPD();
2888 24 : vertices[1] = esd->GetPrimaryVertexTracks();
2889 :
2890 24 : if (!vertices[0]->GetStatus())
2891 : {
2892 0 : AliDebugClass(AliLog::kDebug, "No SPD vertex. Not able to make a reliable multiplicity estimate.");
2893 0 : return -1;
2894 : }
2895 :
2896 24 : if (vertices[0]->IsFromVertexerZ() && vertices[0]->GetDispersion() > 0.02)
2897 : {
2898 0 : AliDebugClass(AliLog::kDebug, "Vertexer z dispersion > 0.02. Not able to make a reliable multiplicity estimate.");
2899 0 : return -2;
2900 : }
2901 :
2902 : Int_t multiplicityEstimate = 0;
2903 :
2904 : // SPD tracklet-only estimate
2905 24 : if (trackType == kTracklets)
2906 : {
2907 0 : const AliMultiplicity* spdmult = esd->GetMultiplicity(); // spd multiplicity object
2908 0 : for (Int_t i=0; i<spdmult->GetNumberOfTracklets(); ++i)
2909 : {
2910 0 : if (TMath::Abs(spdmult->GetEta(i)-etaCent) > etaRange)
2911 : continue; // eta selection for tracklets
2912 0 : multiplicityEstimate++;
2913 0 : }
2914 : return multiplicityEstimate;
2915 : }
2916 :
2917 24 : if (!vertices[1]->GetStatus())
2918 : {
2919 0 : AliDebugClass(AliLog::kDebug, "No track vertex. Not able to make a reliable multiplicity estimate.");
2920 0 : return -3;
2921 : }
2922 :
2923 : // TODO value of displacement to be studied
2924 : const Float_t maxDisplacement = 0.5;
2925 : //check for displaced vertices
2926 24 : Double_t displacement = TMath::Abs(vertices[0]->GetZ() - vertices[1]->GetZ());
2927 24 : if (displacement > maxDisplacement)
2928 : {
2929 0 : AliDebugClass(AliLog::kDebug, Form("Displaced vertices %f > %f",displacement,maxDisplacement));
2930 0 : return -4;
2931 : }
2932 :
2933 : // update eta range in track cuts
2934 24 : float etaMin = etaCent - etaRange, etaMax = etaCent + etaRange;
2935 24 : GetMultEstTrackCuts(kMultEstTrackCutITSSA)->SetEtaRange(etaMin, etaMax);
2936 24 : GetMultEstTrackCuts(kMultEstTrackCutDCAwSPD)->SetEtaRange(etaMin, etaMax);
2937 24 : GetMultEstTrackCuts(kMultEstTrackCutDCAwoSPD)->SetEtaRange(etaMin, etaMax);
2938 :
2939 : //*******************************************************************************************************
2940 : //set counters to initial zeros
2941 : Int_t tracksITSTPC = 0; //number of global tracks for a given event
2942 : Int_t tracksITSSA = 0; //number of ITS standalone tracks for a given event
2943 : Int_t tracksITSTPCSA_complementary = 0; //number of ITS standalone tracks complementary to TPC for a given event
2944 : Int_t trackletsITSTPC_complementary = 0;//number of SPD tracklets complementary to global/ITSSA tracks for a given events
2945 : Int_t trackletsITSSA_complementary = 0; //number of SPD tracklets complementary to ITSSA tracks for a given event
2946 :
2947 24 : const Int_t nESDTracks = esd->GetNumberOfTracks();
2948 :
2949 : // flags for secondary and rejected tracks
2950 : const Int_t kRejBit = BIT(15); // set this bit in global tracks if it is rejected by a cut
2951 : const Int_t kSecBit = BIT(16); // set this bit in global tracks if it is secondary according to a cut
2952 :
2953 900 : for(Int_t itracks=0; itracks < nESDTracks; itracks++) {
2954 426 : esd->GetTrack(itracks)->ResetBit(kSecBit|kRejBit); //reset bits used for flagging secondaries and rejected tracks in case they were changed before this analysis
2955 : }
2956 : const Int_t maxid = nESDTracks; // used to define bool array for check multiple associations of tracklets to one track. array starts at 0.
2957 :
2958 : // bit mask for esd tracks, to check if multiple tracklets are associated to it
2959 24 : TBits globalBits(maxid), pureITSBits(maxid);
2960 : // why labels are used with the data? RS
2961 : // Bool_t globalBits[maxid], pureITSBits[maxid];
2962 : // for(Int_t i=0; i<maxid; i++){ // set all bools to false
2963 : // globalBits[i]=kFALSE;
2964 : // pureITSBits[i]=kFALSE;
2965 : // }
2966 :
2967 : //*******************************************************************************************************
2968 : // get multiplicity from global tracks
2969 900 : for(Int_t itracks = 0; itracks < nESDTracks; itracks++) { // flag the tracks
2970 426 : AliESDtrack* track = esd->GetTrack(itracks);
2971 :
2972 : // if track is a secondary from a V0, flag as a secondary
2973 852 : if (track->IsOn(AliESDtrack::kMultInV0)) {
2974 0 : track->SetBit(kSecBit);
2975 0 : continue;
2976 : }
2977 : /* done via proper DCA cut
2978 : //secondary?
2979 : if (track->IsOn(AliESDtrack::kMultSec)) {
2980 : track->SetBit(kSecBit);
2981 : continue;
2982 : }
2983 : */
2984 : // check tracks with ITS part
2985 : //*******************************************************************************************************
2986 1440 : if (track->IsOn(AliESDtrack::kITSin) && !track->IsOn(AliESDtrack::kITSpureSA) && trackType == kTrackletsITSTPC) { // track has ITS part but is not an ITS_SA
2987 : //*******************************************************************************************************
2988 : // TPC+ITS
2989 588 : if (track->IsOn(AliESDtrack::kTPCin)) { // Global track, has ITS and TPC contributions
2990 492 : if (fgMultEstTrackCuts[kMultEstTrackCutGlobal]->AcceptTrack(track)) { // good ITSTPC track
2991 480 : if (fgMultEstTrackCuts[kMultEstTrackCutDCAwSPD]->AcceptTrack(track) || fgMultEstTrackCuts[kMultEstTrackCutDCAwoSPD]->AcceptTrack(track)) {
2992 222 : tracksITSTPC++; //global track counted
2993 222 : globalBits.SetBitNumber(itracks);
2994 : }
2995 6 : else track->SetBit(kSecBit); // large DCA -> secondary, don't count either track not associated tracklet
2996 : }
2997 18 : else track->SetBit(kRejBit); // bad quality, don't count the track, but may count tracklet if associated
2998 : }
2999 : //*******************************************************************************************************
3000 : // ITS complementary
3001 96 : else if (fgMultEstTrackCuts[kMultEstTrackCutITSSA]->AcceptTrack(track)) { // good ITS complementary track
3002 108 : if (fgMultEstTrackCuts[kMultEstTrackCutDCAwSPD]->AcceptTrack(track) || fgMultEstTrackCuts[kMultEstTrackCutDCAwoSPD]->AcceptTrack(track)) {
3003 48 : tracksITSTPCSA_complementary++;
3004 48 : globalBits.SetBitNumber(itracks);
3005 : }
3006 0 : else track->SetBit(kSecBit); // large DCA -> secondary, don't count either track not associated tracklet
3007 : }
3008 0 : else track->SetBit(kRejBit); // bad quality, don't count the track, but may count tracklet if associated
3009 : }
3010 : //*******************************************************************************************************
3011 : // check tracks from ITS_SA_PURE
3012 1440 : if (track->IsOn(AliESDtrack::kITSin) && track->IsOn(AliESDtrack::kITSpureSA) && trackType == kTrackletsITSSA){
3013 0 : if (fgMultEstTrackCuts[kMultEstTrackCutITSSA]->AcceptTrack(track)) { // good ITSSA track
3014 0 : if (fgMultEstTrackCuts[kMultEstTrackCutDCAwSPD]->AcceptTrack(track) || fgMultEstTrackCuts[kMultEstTrackCutDCAwoSPD]->AcceptTrack(track)) {
3015 0 : tracksITSSA++;
3016 0 : pureITSBits.SetBitNumber(itracks);
3017 : }
3018 0 : else track->SetBit(kSecBit);
3019 : }
3020 0 : else track->SetBit(kRejBit);
3021 : }
3022 426 : }//ESD tracks counted
3023 :
3024 : //*******************************************************************************************************
3025 : // get multiplicity from ITS tracklets to complement TPC+ITS, and ITSpureSA
3026 24 : const AliMultiplicity* spdmult = esd->GetMultiplicity(); // spd multiplicity object
3027 504 : for (Int_t i=0; i<spdmult->GetNumberOfTracklets(); ++i) {
3028 288 : if (TMath::Abs(spdmult->GetEta(i)-etaCent) > etaRange) continue; // eta selection for tracklets
3029 :
3030 : // if counting tracks+tracklets, check if clusters were already used in tracks
3031 144 : Int_t id1, id2, id3, id4;
3032 144 : spdmult->GetTrackletTrackIDs ( i, 0, id1, id2 ); // references for eventual Global/ITS_SA tracks
3033 144 : spdmult->GetTrackletTrackIDs ( i, 1, id3, id4 ); // references for eventual ITS_SA_pure tracks
3034 :
3035 : // are both clusters from the same tracks? If not, skip the tracklet (shouldn't change things much)
3036 288 : if ( ( id1 != id2 && id1 >= 0 && id2 >= 0 ) || ( id3 != id4 && id3 >= 0 && id4 >= 0 ) ) continue;
3037 :
3038 : //referenced track
3039 : //at this point we either have id1 = id2 (id3 = id4) or only one of the ids pair is -1
3040 : // id1>=0, id2>=0, id1=id2 : tracklet has associated track
3041 : // id1>=0, id2 = -1 : 1st layer cluster has associated track
3042 : // id1=-1, id2>=0 : 2nd layer cluster has associated track
3043 : // id1=-1, id2=-1 : tracklet has no associated track
3044 : //
3045 : Int_t bUsedInGlobal(-1);
3046 288 : if ( id1 != -1 ) bUsedInGlobal = globalBits.TestBitNumber(id1) ? id1 : -1;
3047 0 : else if ( id2 != -1) bUsedInGlobal = globalBits.TestBitNumber(id2) ? id2 : -1;
3048 : Int_t bUsedInPureITS(-1);
3049 144 : if ( id3 != -1 ) bUsedInPureITS = pureITSBits.TestBitNumber(id3) ? id3 : -1;
3050 144 : else if ( id4 != -1) bUsedInPureITS = pureITSBits.TestBitNumber(id4) ? id4 : -1;
3051 : //
3052 432 : AliESDtrack* tr_global = bUsedInGlobal >= 0 ? esd->GetTrack ( bUsedInGlobal ) : 0;
3053 288 : AliESDtrack* tr_itssa = bUsedInPureITS >= 0 ? esd->GetTrack ( bUsedInPureITS ) : 0;
3054 : //
3055 : // has associated pure ITS track been associated to a previous tracklet?
3056 : //*******************************************************************************************************
3057 144 : if (trackType == kTrackletsITSTPC) {
3058 : //*******************************************************************************************************
3059 : // count tracklets towards global+complimentary tracks
3060 432 : if ( ( tr_global && !tr_global->TestBit ( kSecBit ) ) && ( tr_global && tr_global->TestBit ( kRejBit ) ) ) { // count tracklet as bad quality track
3061 0 : globalBits.SetBitNumber( bUsedInGlobal ); // mark global track linked to this tracklet as used
3062 0 : ++trackletsITSTPC_complementary;
3063 0 : }
3064 :
3065 144 : if ( bUsedInGlobal < 0 ) ++trackletsITSTPC_complementary; //no associated track, just count the tracklet
3066 : } else {
3067 : //*******************************************************************************************************
3068 : // count tracklets towards ITS_SA_pure tracks
3069 0 : if ( ( tr_itssa && !tr_itssa->TestBit ( kSecBit ) ) && ( tr_itssa && tr_itssa->TestBit ( kRejBit ) ) ) { // count tracklet as bad quality track
3070 0 : pureITSBits.SetBitNumber( bUsedInPureITS ); // mark ITS pure SA track linked to this tracklet as used
3071 0 : ++trackletsITSSA_complementary;
3072 0 : }
3073 :
3074 0 : if ( bUsedInPureITS < 0 ) ++trackletsITSSA_complementary; //no associated track, just count the tracklet
3075 : }
3076 288 : }
3077 :
3078 : //*******************************************************************************************************
3079 24 : if (trackType == kTrackletsITSTPC)
3080 24 : multiplicityEstimate = tracksITSTPC + tracksITSTPCSA_complementary + trackletsITSTPC_complementary;
3081 : else
3082 0 : multiplicityEstimate = tracksITSSA + trackletsITSSA_complementary;
3083 :
3084 : return multiplicityEstimate;
3085 48 : }
3086 :
3087 : //____________________________________________________________________
3088 : void AliESDtrackCuts::SetRequireStandardTOFmatchCuts(){
3089 :
3090 : // setting the TOF cuts flags (kTOFout = TOF matching distance) to true, to include the selection on the standard TOF matching
3091 :
3092 0 : SetRequireTOFout(kTRUE);
3093 0 : SetFlagCutTOFdistance(kTRUE);
3094 0 : SetCutTOFdistance(3.);
3095 :
3096 0 : }
3097 :
|