]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWGPP/TRD/info/AliTRDv0Info.cxx
various extensions
[u/mrichter/AliRoot.git] / PWGPP / TRD / info / AliTRDv0Info.cxx
CommitLineData
1ee39b3a 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: AliTRDv0Info.cxx 27496 2008-07-22 08:35:45Z cblume $ */
17
18////////////////////////////////////////////////////////////////////////////
19// //
20// Reconstruction QA //
21// //
22// Gathers all information necessary for reference data selection about //
23// the track and (in case) its corresponding V0. //
24// Carries out the selection of electrons (from gamma conversions), //
25// pions (from K0s decays) and protons (from Lambda and Anti-Lambda //
26// decays) by cuts specific for the respective decay and particle //
27// species. //
28// (M.Heide, 2009/10/06) //
29// //
30// Authors: //
31// Alex Bercuci <A.Bercuci@gsi.de> //
32// Alex Wilk <wilka@uni-muenster.de> //
33// Markus Heide <mheide@uni-muenster.de> //
34// //
35////////////////////////////////////////////////////////////////////////////
0ed6f095 36
1ee39b3a 37#include "TMath.h"
0ed6f095 38#include "TDatabasePDG.h"
1ee39b3a 39
40#include "AliESDtrack.h"
41#include "AliESDv0.h"
1ee39b3a 42#include "AliLog.h"
0ed6f095 43#include "TVector3.h"
44#include "AliKFParticle.h"
45#include "AliKFVertex.h"
1ee39b3a 46
47#include "AliTRDv0Info.h"
48#include "AliTRDtrackInfo.h"
49#include "AliTRDtrackInfo.h"
50
51ClassImp(AliTRDv0Info)
52
53//_________________________________________________
54AliTRDv0Info::AliTRDv0Info()
55 : TObject()
3d19c1b0 56 ,fQuality(0)
1ee39b3a 57 ,fDCA(10)
58 ,fPointingAngle(10)
59 ,fOpenAngle(10)
60 ,fPsiPair(99)
61 ,fMagField(0)
62 ,fRadius(0)
1ee39b3a 63 ,fV0Momentum(0)
1ee39b3a 64 ,fNindex(0)
65 ,fPindex(0)
0ed6f095 66 ,fInputEvent(NULL)
67 ,fPrimaryVertex(NULL)
64d57299 68 ,fTrackP(NULL)
69 ,fTrackN(NULL)
1ee39b3a 70{
71 //
72 // Default constructor
73 //
74
1ee39b3a 75 memset(fDetPID, 0, 2*kNDaughters*kNDetectors*AliPID::kSPECIES*sizeof(Float_t));
8bc8ea55 76 memset(fComPID, 0, 2*kNDaughters*AliPID::kSPECIES*sizeof(Float_t));
61cfa442 77 memset(fInvMass, 0, kNDecays*sizeof(Double_t));
0ed6f095 78 memset(fArmenteros, 0, kNDecays*sizeof(Bool_t));
61cfa442 79 memset(fTPCdEdx, 0, kNDaughters*sizeof(Float_t));
0ed6f095 80 memset(fChi2ndf, 0, kNDecays*sizeof(Double_t));
7fe4e88b 81 memset(fDownOpenAngle, 0, kNDecays*sizeof(Float_t));
82 memset(fDownPsiPair, 0, kNDecays*sizeof(Float_t));
1ee39b3a 83 /////////////////////////////////////////////////////////////////////////////
84 //Set Cut values: First specify decay in brackets, then the actual cut value!
85 /////////////////////////////////////////////////////////////////////////////
86
87 //Upper limit for distance of closest approach of two daughter tracks :
8bc8ea55 88 fUpDCA[kGamma] = 1000.;
89 fUpDCA[kK0s] = 0.08;
90 fUpDCA[kLambda] = 0.2;
91 fUpDCA[kAntiLambda] = 0.2;
1ee39b3a 92
93 //Upper limit for pointing angle (= angle between between vector from primary to secondary vertex and reconstructed momentum of V0 mother particle) :
94 fUpPointingAngle[kGamma] = 0.03;
95 fUpPointingAngle[kK0s] = 0.03;
8bc8ea55 96 fUpPointingAngle[kLambda] = 0.04;
97 fUpPointingAngle[kAntiLambda] = 0.04;
1ee39b3a 98
99 //Upper limit for invariant mass of V0 mother :
8bc8ea55 100 fUpInvMass[kGamma][0] = 0.05;// second pair of brackets is for momentum bin: 0: below mother momentm of 2.5 GeV
1ee39b3a 101 fUpInvMass[kGamma][1] = 0.07;//1: above 2.5 GeV
8bc8ea55 102 fUpInvMass[kK0s][0] = fUpInvMass[kK0s][1] = 0.50265;
103 fUpInvMass[kLambda][0] = fUpInvMass[kLambda][1] = 1.1207;
104 fUpInvMass[kAntiLambda][0] = fUpInvMass[kAntiLambda][1] = 1.1207;
1ee39b3a 105
106 //Lower limit for invariant mass of V0 mother :
107 fDownInvMass[kGamma] = -1.;
8bc8ea55 108 fDownInvMass[kK0s] = 0.49265;
109 fDownInvMass[kLambda] = 1.107;
110 fDownInvMass[kAntiLambda] = 1.107;
1ee39b3a 111
0ed6f095 112 //Upper limit for KF Chi2/NDF value;
113 fUpChi2ndf[kGamma] = 10000.;//7.;
114 fUpChi2ndf[kK0s] = 10000.;//5.;
115 fUpChi2ndf[kLambda] = 10000.;//5.;
116 fUpChi2ndf[kAntiLambda] = 10000.;//5.;
117
1ee39b3a 118 //Lower limit for distance from secondary vertex to primary vertex in x-y plane :
8bc8ea55 119 fDownRadius[kGamma] = 6.;
1ee39b3a 120 fDownRadius[kK0s] = 0.;
8bc8ea55 121 fDownRadius[kLambda] = 0.;
122 fDownRadius[kAntiLambda] = 0.;
1ee39b3a 123
124 //Upper limit for distance from secondary vertex to primary vertex in x-y plane :
125 fUpRadius[kGamma] = 1000.;
8bc8ea55 126 fUpRadius[kK0s] = 20.;
1ee39b3a 127 fUpRadius[kLambda] = 1000.;
128 fUpRadius[kAntiLambda] = 1000.;
129
130 //Upper limit for opening angle between two daughter tracks (characteristically near zero for conversions) :
131 fUpOpenAngle[kGamma] = 0.1;
132 fUpOpenAngle[kK0s] = 3.15;
133 fUpOpenAngle[kLambda] = 3.15;
134 fUpOpenAngle[kAntiLambda] = 3.15;
135
136 //Upper limit for angle between daughter momentum plane and plane perpendicular to magnetic field (characteristically around zero for conversions) :
8bc8ea55 137 fUpPsiPair[kGamma] = 0.05;
1ee39b3a 138 fUpPsiPair[kK0s] = 1.6;
139 fUpPsiPair[kLambda] = 1.6;
140 fUpPsiPair[kAntiLambda] = 1.6;
141
142 //Lower limit for likelihood value of TPC PID :
0ed6f095 143 fDownTPCPIDneg[AliPID::kElectron] = 0.;
144 fDownTPCPIDpos[AliPID::kElectron] = 0.;
145
146 fDownTPCPIDneg[AliPID::kMuon] = 0.;
147 fDownTPCPIDpos[AliPID::kMuon] = 0.;
148
149 fDownTPCPIDneg[AliPID::kPion] = 0.;
150 fDownTPCPIDpos[AliPID::kPion] = 0.;
151
152 fDownTPCPIDneg[AliPID::kKaon] = 0.;
153 fDownTPCPIDpos[AliPID::kKaon] = 0.;
154
155 fDownTPCPIDneg[AliPID::kProton] = 0.;
156 fDownTPCPIDpos[AliPID::kProton] = 0.;
157
158 //Lower limit for likelihood value of combined PID :
159 fDownComPIDneg[AliPID::kElectron] = 0.;
160 fDownComPIDpos[AliPID::kElectron] = 0.;
161
162 fDownComPIDneg[AliPID::kMuon] = 0.;
163 fDownComPIDpos[AliPID::kMuon] = 0.;
164
165 fDownComPIDneg[AliPID::kPion] = 0.;
166 fDownComPIDpos[AliPID::kPion] = 0.;
167
168 fDownComPIDneg[AliPID::kKaon] = 0.;
169 fDownComPIDpos[AliPID::kKaon] = 0.;
170
171 fDownComPIDneg[AliPID::kProton] = 0.;
172 fDownComPIDpos[AliPID::kProton] = 0.;
173
174 //Lower limit for likelihood value of combined PID for daughter track which doesn't enter reference data (here: pion daughters from Lambda decays:
175 fDownComPIDnegPart[AliPID::kElectron] = 0.;
176 fDownComPIDposPart[AliPID::kElectron] = 0.;
177
178 fDownComPIDnegPart[AliPID::kMuon] = 0.;
179 fDownComPIDposPart[AliPID::kMuon] = 0.;
180
181 fDownComPIDnegPart[AliPID::kPion] = 0.;
182 fDownComPIDposPart[AliPID::kPion] = 0.;
183
184 fDownComPIDnegPart[AliPID::kKaon] = 0.;
185 fDownComPIDposPart[AliPID::kKaon] = 0.;
186
187 fDownComPIDnegPart[AliPID::kProton] = 0.;
188 fDownComPIDposPart[AliPID::kProton] = 0.;
189
190 //Parameters for data with well-calibrated PID (after usage of tender):
191 /* //Lower limit for likelihood value of TPC PID :
1ee39b3a 192 fDownTPCPIDneg[AliPID::kElectron] = 0.21;
193 fDownTPCPIDpos[AliPID::kElectron] = 0.21;
194
195 fDownTPCPIDneg[AliPID::kMuon] = 0.21;
196 fDownTPCPIDpos[AliPID::kMuon] = 0.21;
197
198 fDownTPCPIDneg[AliPID::kPion] = 0.21;
199 fDownTPCPIDpos[AliPID::kPion] = 0.21;
200
201 fDownTPCPIDneg[AliPID::kKaon] = 0.21;
202 fDownTPCPIDpos[AliPID::kKaon] = 0.21;
203
204 fDownTPCPIDneg[AliPID::kProton] = 0.21;
205 fDownTPCPIDpos[AliPID::kProton] = 0.21;
8bc8ea55 206
0ed6f095 207 //Lower limit for likelihood value of combined PID :
8bc8ea55 208 fDownComPIDneg[AliPID::kElectron] = 0.21;
209 fDownComPIDpos[AliPID::kElectron] = 0.21;
210
211 fDownComPIDneg[AliPID::kMuon] = 0.21;
212 fDownComPIDpos[AliPID::kMuon] = 0.21;
213
214 fDownComPIDneg[AliPID::kPion] = 0.9;
215 fDownComPIDpos[AliPID::kPion] = 0.9;
216
217 fDownComPIDneg[AliPID::kKaon] = 0.21;
218 fDownComPIDpos[AliPID::kKaon] = 0.21;
219
220 fDownComPIDneg[AliPID::kProton] = 0.9;
221 fDownComPIDpos[AliPID::kProton] = 0.9;
222
223 //Lower limit for likelihood value of combined PID for daughter track which doesn't enter reference data (here: pion daughters from Lambda decays:
224 fDownComPIDnegPart[AliPID::kElectron] = 0.05;
225 fDownComPIDposPart[AliPID::kElectron] = 0.05;
226
227 fDownComPIDnegPart[AliPID::kMuon] = 0.05;
228 fDownComPIDposPart[AliPID::kMuon] = 0.05;
229
230 fDownComPIDnegPart[AliPID::kPion] = 0.05;
231 fDownComPIDposPart[AliPID::kPion] = 0.05;
232
233 fDownComPIDnegPart[AliPID::kKaon] = 0.05;
234 fDownComPIDposPart[AliPID::kKaon] = 0.05;
235
236 fDownComPIDnegPart[AliPID::kProton] = 0.05;
0ed6f095 237 fDownComPIDposPart[AliPID::kProton] = 0.05;*/
3d19c1b0 238}
239
240//_________________________________________________
241AliTRDv0Info::AliTRDv0Info(const AliTRDv0Info &ref)
5047978d 242 : TObject((TObject&)ref)
3d19c1b0 243 ,fQuality(ref.fQuality)
244 ,fDCA(ref.fDCA)
245 ,fPointingAngle(ref.fPointingAngle)
246 ,fOpenAngle(ref.fOpenAngle)
247 ,fPsiPair(ref.fPsiPair)
248 ,fMagField(ref.fMagField)
0ed6f095 249 ,fRadius(ref.fRadius)
3d19c1b0 250 ,fV0Momentum(ref.fV0Momentum)
3d19c1b0 251 ,fNindex(ref.fNindex)
252 ,fPindex(ref.fPindex)
0ed6f095 253 ,fInputEvent(ref.fInputEvent)
254 ,fPrimaryVertex(ref.fPrimaryVertex)
64d57299 255 ,fTrackP(ref.fTrackP)
256 ,fTrackN(ref.fTrackN)
3d19c1b0 257{
258 //
259 // Copy constructor
260 //
64d57299 261
3d19c1b0 262 memcpy(fDetPID, ref.fDetPID, 2*kNDaughters*kNDetectors*AliPID::kSPECIES*sizeof(Float_t));
263 memcpy(fComPID, ref.fComPID, 2*kNDaughters*AliPID::kSPECIES*sizeof(Float_t));
61cfa442 264 memcpy(fInvMass, ref.fInvMass, kNDecays*sizeof(Double_t));
0ed6f095 265 memcpy(fArmenteros, ref.fArmenteros, kNDecays*sizeof(Bool_t));
266 memcpy(fChi2ndf, ref.fChi2ndf, kNDecays*sizeof(Double_t));
267 memcpy(fTPCdEdx, ref.fTPCdEdx, kNDaughters*sizeof(Float_t));
1ee39b3a 268
3d19c1b0 269 //Upper limit for distance of closest approach of two daughter tracks :
270 memcpy(fUpDCA, ref.fUpDCA, kNDecays*sizeof(Float_t));
271 memcpy(fUpPointingAngle, ref.fUpPointingAngle, kNDecays*sizeof(Float_t));
272 memcpy(fUpOpenAngle, ref.fUpOpenAngle, kNDecays*sizeof(Float_t));
273 memcpy(fDownOpenAngle, ref.fDownOpenAngle, kNDecays*sizeof(Float_t));
274 memcpy(fUpPsiPair, ref.fUpPsiPair, kNDecays*sizeof(Float_t));
275 memcpy(fDownPsiPair, ref.fDownPsiPair, kNDecays*sizeof(Float_t));
276 memcpy(fUpInvMass, ref.fUpInvMass, kNDecays*kNMomBins*sizeof(Double_t));
277 memcpy(fDownInvMass, ref.fDownInvMass, kNDecays*sizeof(Double_t));
0ed6f095 278 memcpy(fUpChi2ndf, ref.fUpChi2ndf, kNDecays*sizeof(Double_t));
3d19c1b0 279 memcpy(fUpRadius, ref.fUpRadius, kNDecays*sizeof(Float_t));
280 memcpy(fDownRadius, ref.fDownRadius, kNDecays*sizeof(Float_t));
281 memcpy(fDownTPCPIDneg, ref.fDownTPCPIDneg, AliPID::kSPECIES*sizeof(Float_t));
282 memcpy(fDownTPCPIDpos, ref.fDownTPCPIDpos, AliPID::kSPECIES*sizeof(Float_t));
283 memcpy(fDownComPIDneg, ref.fDownComPIDneg, AliPID::kSPECIES*sizeof(Float_t));
284 memcpy(fDownComPIDpos, ref.fDownComPIDpos, AliPID::kSPECIES*sizeof(Float_t));
285 memcpy(fDownComPIDnegPart, ref.fDownComPIDnegPart, AliPID::kSPECIES*sizeof(Float_t));
286 memcpy(fDownComPIDposPart, ref.fDownComPIDposPart, AliPID::kSPECIES*sizeof(Float_t));
1ee39b3a 287}
288
289//_________________________________________________
b37d601d 290void AliTRDv0Info::SetV0Info(const AliESDv0 *esdv0)
64d57299 291{
292 //Gets values of ESDv0 and daughter track properties
1ee39b3a 293 //See header file for description of variables
294
1ee39b3a 295 fQuality = Quality(esdv0);//Attributes an Int_t to the V0 due to quality cuts (= 1 if V0 is accepted, other integers depending on cut which excludes the vertex)
296
297 fRadius = Radius(esdv0);//distance from secondary vertex to primary vertex in x-y plane
298
299 fDCA = esdv0->GetDcaV0Daughters();//distance of closest approach of two daughter tracks
300
301 fPointingAngle = TMath::ACos(esdv0->GetV0CosineOfPointingAngle());// pointing angle (= angle between between vector from primary to secondary vertex and reconstructed momentum of V0 mother particle)
302
303 fOpenAngle = OpenAngle(esdv0);//Opening angle between two daughter tracks
304
305 fPsiPair = PsiPair(esdv0);//Angle between daughter momentum plane and plane perpendicular to magnetic field
306
307 fV0Momentum = V0Momentum(esdv0);//Reconstructed momentum of the mother particle
308
3d19c1b0 309 //4 decay types : conversions, K0s, Lambda, Anti-Lambda
1ee39b3a 310 //five particle types: electrons, muons, pions, kaons, protons (muons and kaons not involved)
3d19c1b0 311 for(Int_t idecay(0), part1(-1), part2(-1); idecay < kNDecays; idecay++){
0ed6f095 312
313 fArmenteros[idecay]=Armenteros(esdv0, idecay);//Attribute the Armenteros yes/no decision for every decay type
a120f5e2 314 switch(idecay){
315 case kLambda: //protons and pions from Lambda
3d19c1b0 316 part1 = AliPID::kProton;
317 part2 = AliPID::kPion;
a120f5e2 318 break;
319 case kAntiLambda: //antiprotons and pions from Anti-Lambda
3d19c1b0 320 part1 = AliPID::kPion;
321 part2 = AliPID::kProton;
a120f5e2 322 break;
323 case kK0s: //pions from K0s
3d19c1b0 324 part1 = part2 = AliPID::kPion;
a120f5e2 325 break;
326 case kGamma://electrons from conversions
3d19c1b0 327 part1 = part2 = AliPID::kElectron;
a120f5e2 328 break;
329 }
3d19c1b0 330 fInvMass[idecay] = InvMass(part1, part2, esdv0);//Calculate invariant mass for all of our four supposed decays
24b3cfb9 331
332 // Comment out until bug fix is provided
333 // A.Bercuci 14. July 2010
334 //fChi2ndf[idecay] = KFChi2ndf(part1, part2,idecay);
0ed6f095 335
1ee39b3a 336 }
3d19c1b0 337 //Gets all likelihood values from TPC, TOF and ITS PID for the fDetPID[kNDaughters][kNDetectors][AliPID::kSPECIES] array
338 GetDetectorPID();
339 //Bayesian combination of likelihoods from TPC and TOF
340 CombinePID();
0ed6f095 341 //TPC dE/dx values for both tracks
342 GetTPCdEdx();
343
1ee39b3a 344}
345//_________________________________________________
b37d601d 346Float_t AliTRDv0Info::V0Momentum(const AliESDv0 *esdv0) const
1ee39b3a 347{
348 //
349 // Reconstructed momentum of V0 mother particle
350 //
351
352 Double_t mn[3] = {0,0,0};
353 Double_t mp[3] = {0,0,0};
354
355
356 esdv0->GetNPxPyPz(mn[0],mn[1],mn[2]);//reconstructed cartesian momentum components of negative daughter;
357 esdv0->GetPPxPyPz(mp[0],mp[1],mp[2]);//reconstructed cartesian momentum components of positive daughter;
358
359
360 return TMath::Sqrt((mn[0]+mp[0])*(mn[0]+mp[0]) + (mn[1]+mp[1])*(mn[1]+mp[1])+(mn[2]+mp[2])*(mn[2]+mp[2]));
361}
362
363//_________________________________________________
b37d601d 364Double_t AliTRDv0Info::InvMass(Int_t part1, Int_t part2, const AliESDv0 *esdv0) const
1ee39b3a 365{
366 //
367 // Invariant mass of reconstructed V0 mother
368 //
369
370 const Double_t kpmass[5] = {AliPID::ParticleMass(AliPID::kElectron),AliPID::ParticleMass(AliPID::kMuon),AliPID::ParticleMass(AliPID::kPion),AliPID::ParticleMass(AliPID::kKaon),AliPID::ParticleMass(AliPID::kProton)};
371 //Masses of electrons, muons, pions, kaons and protons, as implemented in ROOT
372
373
374 Double_t mn[3] = {0,0,0};
375 Double_t mp[3] = {0,0,0};
376
377 esdv0->GetNPxPyPz(mn[0],mn[1],mn[2]);//reconstructed cartesian momentum components of negative daughter;
378 esdv0->GetPPxPyPz(mp[0],mp[1],mp[2]);//reconstructed cartesian momentum components of positive daughter;
379
0ed6f095 380 Double_t mass1 = kpmass[part1];//sets supposed rest masses for both daughters: positive
381 Double_t mass2 = kpmass[part2];//negative
1ee39b3a 382
383 //Calculate daughters' energies :
384 Double_t e1 = TMath::Sqrt(mass1*mass1+
385 mp[0]*mp[0]+
386 mp[1]*mp[1]+
387 mp[2]*mp[2]);
388 Double_t e2 = TMath::Sqrt(mass2*mass2+
389 mn[0]*mn[0]+
390 mn[1]*mn[1]+
391 mn[2]*mn[2]);
392
393 //Sum of daughter momenta :
394 Double_t momsum =
395 (mn[0]+mp[0])*(mn[0]+mp[0])+
396 (mn[1]+mp[1])*(mn[1]+mp[1])+
397 (mn[2]+mp[2])*(mn[2]+mp[2]);
398
399 //invariant mass :
400 Double_t mInv = TMath::Sqrt((e1+e2)*(e1+e2)-momsum);
401
402 return mInv;
403
404}
405//_________________________________________________
b37d601d 406Float_t AliTRDv0Info::OpenAngle(const AliESDv0 *esdv0)
64d57299 407{
408 //Opening angle between two daughter tracks
1ee39b3a 409 Double_t mn[3] = {0,0,0};
410 Double_t mp[3] = {0,0,0};
411
412
413 esdv0->GetNPxPyPz(mn[0],mn[1],mn[2]);//reconstructed cartesian momentum components of negative daughter;
414 esdv0->GetPPxPyPz(mp[0],mp[1],mp[2]);//reconstructed cartesian momentum components of positive daughter;
415
416
417 fOpenAngle = TMath::ACos((mp[0]*mn[0] + mp[1]*mn[1] + mp[2]*mn[2])/(TMath::Sqrt(mp[0]*mp[0] + mp[1]*mp[1] + mp[2]*mp[2])*TMath::Sqrt(mn[0]*mn[0] + mn[1]*mn[1] + mn[2]*mn[2])));
418
419 return fOpenAngle;
420}
421
422//_________________________________________________
b37d601d 423Float_t AliTRDv0Info::PsiPair(const AliESDv0 *esdv0)
64d57299 424{
425 //Angle between daughter momentum plane and plane perpendicular to magnetic field
1ee39b3a 426 Double_t x, y, z;
427 esdv0->GetXYZ(x,y,z);//Reconstructed coordinates of V0; to be replaced by Markus Rammler's method in case of conversions!
428
429 Double_t mn[3] = {0,0,0};
430 Double_t mp[3] = {0,0,0};
431
432
433 esdv0->GetNPxPyPz(mn[0],mn[1],mn[2]);//reconstructed cartesian momentum components of negative daughter;
434 esdv0->GetPPxPyPz(mp[0],mp[1],mp[2]);//reconstructed cartesian momentum components of positive daughter;
435
436
437 Double_t deltat = 1.;
438 deltat = TMath::ATan(mp[2]/(TMath::Sqrt(mp[0]*mp[0] + mp[1]*mp[1])+1.e-13)) - TMath::ATan(mn[2]/(TMath::Sqrt(mn[0]*mn[0] + mn[1]*mn[1])+1.e-13));//difference of angles of the two daughter tracks with z-axis
439
440 Double_t radiussum = TMath::Sqrt(x*x + y*y) + 50;//radius to which tracks shall be propagated
441
442 Double_t momPosProp[3];
443 Double_t momNegProp[3];
444
445 AliExternalTrackParam nt(*fTrackN), pt(*fTrackP);
446
447 fPsiPair = 4.;
448
449 if(nt.PropagateTo(radiussum,fMagField) == 0)//propagate tracks to the outside
450 fPsiPair = -5.;
451 if(pt.PropagateTo(radiussum,fMagField) == 0)
452 fPsiPair = -5.;
453 pt.GetPxPyPz(momPosProp);//Get momentum vectors of tracks after propagation
454 nt.GetPxPyPz(momNegProp);
455
456 Double_t pEle =
457 TMath::Sqrt(momNegProp[0]*momNegProp[0]+momNegProp[1]*momNegProp[1]+momNegProp[2]*momNegProp[2]);//absolute momentum value of negative daughter
458 Double_t pPos =
459 TMath::Sqrt(momPosProp[0]*momPosProp[0]+momPosProp[1]*momPosProp[1]+momPosProp[2]*momPosProp[2]);//absolute momentum value of positive daughter
460
461 Double_t scalarproduct =
462 momPosProp[0]*momNegProp[0]+momPosProp[1]*momNegProp[1]+momPosProp[2]*momNegProp[2];//scalar product of propagated positive and negative daughters' momenta
463
464 Double_t chipair = TMath::ACos(scalarproduct/(pEle*pPos));//Angle between propagated daughter tracks
465
466 fPsiPair = TMath::Abs(TMath::ASin(deltat/chipair));
467
468 return fPsiPair;
469
470}
0ed6f095 471//_________________________________________________
472Double_t AliTRDv0Info::KFChi2ndf(Int_t part1, Int_t part2,Int_t decay){
473 //Calculates Kalman filter Chi2/NDF
474 Int_t mothers[4]={22,310,3122,3122};
475
476 const Double_t partMass=TDatabasePDG::Instance()->GetParticle(mothers[decay])->Mass();
477 const Double_t massWidth[4] = {0.001, 0., 0., 0.};
478
479 AliKFParticle *kfMother = CreateMotherParticle(fTrackP, fTrackN, part1, part2);
480
481 // Lambda
482 if(!kfMother) {
483 return kFALSE;
484 }
485
486 // production vertex is set in the 'CreateMotherParticle' function
487 kfMother->SetMassConstraint(partMass, massWidth[decay]);
488
489 Double_t chi2ndf = (kfMother->GetChi2()/kfMother->GetNDF());
490
491 if(kfMother)delete kfMother;
492 return chi2ndf;
493}
494//________________________________________________________________
64d57299 495AliKFParticle *AliTRDv0Info::CreateMotherParticle(const AliESDtrack *pdaughter, const AliESDtrack *ndaughter, Int_t pspec, Int_t nspec){
0ed6f095 496 //
2c33fb46 497 // Creates a mother particle on the HEAP !! User code is responsible for its deletion
0ed6f095 498 //
499 AliKFParticle pkfdaughter(*pdaughter, pspec);
500 AliKFParticle nkfdaughter(*ndaughter, nspec);
501
502
503 // Create the mother particle
504 AliKFParticle *m = new AliKFParticle(pkfdaughter, nkfdaughter);
505
506 AliKFVertex improvedVertex = *fPrimaryVertex;
507 improvedVertex += *m;
508 m->SetProductionVertex(improvedVertex);
509
1ee39b3a 510
0ed6f095 511 return m;
512}
1ee39b3a 513//_________________________________________________
64d57299 514Int_t AliTRDv0Info::HasTrack(const AliTRDtrackInfo * const track) const
3d19c1b0 515{
64d57299 516 //Checks if track is a secondary vertex daughter (due to V0 finder)
1ee39b3a 517
d80a6a00 518 if(!track) return 0;
519 if(!fTrackP->GetID()) return 0;
520 if(!fTrackN->GetID()) return 0;
521
3d19c1b0 522 Int_t trackID(track->GetTrackId());//index of the track
b9ddd472 523 return HasTrack(trackID);
524}
1ee39b3a 525
b9ddd472 526//_________________________________________________
64d57299 527Int_t AliTRDv0Info::HasTrack(Int_t trackID) const
b9ddd472 528{
1ee39b3a 529 //comparing index of track with indices of pos./neg. V0 daughter :
b9ddd472 530 if(fNindex==trackID) return -1;
531 else if(fPindex==trackID) return 1;
532 else return 0;
1ee39b3a 533}
3d19c1b0 534
1ee39b3a 535//_________________________________________________
536void AliTRDv0Info::GetDetectorPID()
64d57299 537{
538 //PID likelihoods from TPC, TOF, and ITS, for all particle species
1ee39b3a 539
540 fTrackN->GetTPCpid(fDetPID[kNeg][kTPC]);
541 fTrackP->GetTPCpid(fDetPID[kPos][kTPC]);
542 fTrackN->GetTOFpid(fDetPID[kNeg][kTOF]);
543 fTrackP->GetTOFpid(fDetPID[kPos][kTOF]);
544 fTrackN->GetITSpid(fDetPID[kNeg][kITS]);
545 fTrackP->GetITSpid(fDetPID[kPos][kITS]);
546
8bc8ea55 547 Long_t statusN = fTrackN->GetStatus();
548 Long_t statusP = fTrackP->GetStatus();
549
550 if(!(statusN & AliESDtrack::kTPCpid)){
551 for(Int_t iPart = 0; iPart < AliPID::kSPECIES; iPart++){
552 fDetPID[kNeg][kTPC][iPart] = 0.2;
553 }
554 }
555 if(!(statusN & AliESDtrack::kTOFpid)){
556 for(Int_t iPart = 0; iPart < AliPID::kSPECIES; iPart++){
557 fDetPID[kNeg][kTOF][iPart] = 0.2;
558 }
559
560 }
561 if(!(statusN & AliESDtrack::kITSpid)){
562 for(Int_t iPart = 0; iPart < AliPID::kSPECIES; iPart++){
563 fDetPID[kNeg][kITS][iPart] = 0.2;
564 }
565 }
566 if(!(statusP & AliESDtrack::kTPCpid)){
567 for(Int_t iPart = 0; iPart < AliPID::kSPECIES; iPart++){
568 fDetPID[kPos][kTPC][iPart] = 0.2;
569 }
570 }
571 if(!(statusP & AliESDtrack::kTOFpid)){
572 for(Int_t iPart = 0; iPart < AliPID::kSPECIES; iPart++){
573 fDetPID[kPos][kTOF][iPart] = 0.2;
574 }
575
576 }
577 if(!(statusP & AliESDtrack::kITSpid)){
578 for(Int_t iPart = 0; iPart < AliPID::kSPECIES; iPart++){
579 fDetPID[kPos][kITS][iPart] = 0.2;
580 }
581 }
1ee39b3a 582
8bc8ea55 583}
584//____________________________________________________________________________________
585void AliTRDv0Info::CombinePID()
586{
64d57299 587 //combined bayesian PID from TPC and TOF
8bc8ea55 588 Double_t partrat[AliPID::kSPECIES] = {0.208, 0.010, 0.662, 0.019, 0.101};
589
590 for(Int_t iSign = 0; iSign < kNDaughters; iSign++)
591 {
592 for(Int_t iPart = 0; iPart < AliPID::kSPECIES; iPart++)
593 {
594 fComPID[iSign][iPart] = (partrat[iPart]*fDetPID[iSign][kTPC][iPart]*fDetPID[iSign][kTOF][iPart])/((partrat[0]*fDetPID[iSign][kTPC][0]*fDetPID[iSign][kTOF][0])+(partrat[1]*fDetPID[iSign][kTPC][1]*fDetPID[iSign][kTOF][1])+(partrat[2]*fDetPID[iSign][kTPC][2]*fDetPID[iSign][kTOF][2])+(partrat[3]*fDetPID[iSign][kTPC][3]*fDetPID[iSign][kTOF][3])+(partrat[4]*fDetPID[iSign][kTPC][4]*fDetPID[iSign][kTOF][4]));
595
596 }
597 }
598}
1ee39b3a 599//_________________________________________________
d80a6a00 600Bool_t AliTRDv0Info::GetTPCdEdx()
0ed6f095 601{
64d57299 602 //gets the TPC dE/dx for both daughter tracks
d80a6a00 603 if(!fTrackP->GetID()) return 0;
604 if(!fTrackN->GetID()) return 0;
605
0ed6f095 606 fTPCdEdx[kNeg] = fTrackN->GetTPCsignal();
607 fTPCdEdx[kPos] = fTrackP->GetTPCsignal();
d80a6a00 608 return 1;
0ed6f095 609
610}
611//_________________________________________________
b37d601d 612Bool_t AliTRDv0Info::TPCdEdxCuts(Int_t part, const AliTRDtrackInfo * const track)
0ed6f095 613{
64d57299 614 //applies cuts on TPC dE/dx according to particle species; cutting lines are drawn shifted to the Bethe-Bloch paremeterization
d80a6a00 615 if(!fTrackP->GetID()) return 0;
616 if(!fTrackN->GetID()) return 0;
617
0ed6f095 618 //Bethe-Bloch lines
619 Double_t alephParameters[5];
620
621 // data
622 alephParameters[0] = 0.0283086;
623 alephParameters[1] = 2.63394e+01;
624 alephParameters[2] = 5.04114e-11;
625 alephParameters[3] = 2.12543e+00;
626 alephParameters[4] = 4.88663e+00;
627
628
629 Double_t deposit = 0;
630 Float_t x = 0;
631 if(HasTrack(track) == 1){
632 x = fTrackP->P();
633 deposit = fTPCdEdx[kPos];
634 }
635 else if(HasTrack(track) == -1){
636 x = fTrackN->P();
637 deposit = fTPCdEdx[kNeg];
638 }
639 else{
640 printf("No track found");
641 return 0;
642 }
643 if(x < 0.2)return 0;
644
3e628d57 645 Double_t upLimits[5]={
646 85.,
647 1000.,
648 50.*AliExternalTrackParam::BetheBlochAleph(x/0.13957, alephParameters[0], alephParameters[1], alephParameters[2], alephParameters[3], alephParameters[4])+6.,
649 1000.,
650 50.*AliExternalTrackParam::BetheBlochAleph(x/0.93827, alephParameters[0], alephParameters[1], alephParameters[2], alephParameters[3], alephParameters[4])+10.};
651 Double_t downLimits[5]={
652 62.,
653 40.,
654 50.*AliExternalTrackParam::BetheBlochAleph(x/0.13957, alephParameters[0], alephParameters[1], alephParameters[2], alephParameters[3], alephParameters[4])-6.,
655 40.,
656 50.*AliExternalTrackParam::BetheBlochAleph(x/0.93827, alephParameters[0], alephParameters[1], alephParameters[2], alephParameters[3], alephParameters[4])-11.};
0ed6f095 657
658
659 if(x < 0.7){
660 downLimits[4]=90;
661 }
662 if(x < 1.25){
663 upLimits[0] = 85;
664 }
665 else{
666 downLimits[0] = 64;
667 }
668
669
670 if(deposit < downLimits[part])
671 return 0;
672 if(deposit > upLimits[part])
673 return 0;
674
675
676 return 1;
677
678}
679//_________________________________________________
b37d601d 680Float_t AliTRDv0Info::Radius(const AliESDv0 *esdv0)
64d57299 681{
682 //distance from secondary vertex to primary vertex in x-y plane
1ee39b3a 683 Double_t x, y, z;
0ed6f095 684 esdv0->GetXYZ(x,y,z); //Reconstructed coordinates of V0
1ee39b3a 685 fRadius = TMath::Sqrt(x*x + y*y);
686 return fRadius;
687
688}
689
690//_________________________________________________
b37d601d 691Int_t AliTRDv0Info::Quality(const AliESDv0 *const esdv0)
1ee39b3a 692{
693 //
694 // Checking track and V0 quality status in order to exclude vertices based on poor information
695 //
696
697 Float_t nClsN;
698 nClsN = fTrackN->GetTPCNcls();//number of found clusters in TPC for negative track
699 Float_t nClsFN;
700 nClsFN = fTrackN->GetTPCNclsF();//number of findable clusters in TPC for negative track
701 Float_t nClsP;
702 nClsP = fTrackP->GetTPCNcls();//number of found clusters in TPC for positive track
703 Float_t nClsFP;
704 nClsFP = fTrackP->GetTPCNclsF();//number of findable clusters in TPC for positive track
705
706 fQuality = 0;
707
708
0ed6f095 709 if (!(esdv0->GetOnFlyStatus()))//accept only vertices from online V0 finder
710 return -1;
711
1ee39b3a 712 Float_t clsRatioN;
713 Float_t clsRatioP;
714
0ed6f095 715 if((nClsFN < 80) || (nClsFP < 80)) return -2;//reject all V0s where at least one track has less than 80 TPC clusters
716
717 // Chi2 per TPC cluster
718 Int_t nTPCclustersP = fTrackP->GetTPCclusters(0);
719 Int_t nTPCclustersN = fTrackN->GetTPCclusters(0);
720 Float_t chi2perTPCclusterP = fTrackP->GetTPCchi2()/Float_t(nTPCclustersP);
721 Float_t chi2perTPCclusterN = fTrackN->GetTPCchi2()/Float_t(nTPCclustersN);
722
723 if((chi2perTPCclusterN > 3.5)||(chi2perTPCclusterP > 3.5)) return -3;//reject all V0s where at least one track has a chi2 above 3.5
1ee39b3a 724
725 clsRatioN = nClsN/nClsFN; //ratios of found to findable clusters in TPC
726 clsRatioP = nClsP/nClsFP;
0ed6f095 727
728 if((clsRatioN < 0.6)||(clsRatioP < 0.6))//exclude tracks with low ratio of found to findable TPC clusters
729 return -4;
730
1ee39b3a 731 if (!((fTrackP->GetStatus() &
732 AliESDtrack::kTPCrefit)))//accept only vertices in which both tracks have TPC refit
0ed6f095 733 return -5;
1ee39b3a 734 if (!((fTrackN->GetStatus() &
735 AliESDtrack::kTPCrefit)))
0ed6f095 736 return -6;
1ee39b3a 737 if (fTrackP->GetKinkIndex(0)>0 ||
738 fTrackN->GetKinkIndex(0)>0 )//exclude tracks with kinks
e148d6e6 739 return -7;
0ed6f095 740
741 if(!(V0SignCheck()))
742 return -8;
1ee39b3a 743 fQuality = 1;
744 return fQuality;
745}
5047978d 746
0ed6f095 747//________________________________________________________________
748Bool_t AliTRDv0Info::V0SignCheck(){
749 //
750 // Check if v0 daughters really carry opposite charges
751 //
752
753 Int_t qP = fTrackP->Charge();
754 Int_t qN = fTrackN->Charge();
755
756 if((qP*qN) != -1) return kFALSE;
757
758 return kTRUE;
759}
5047978d 760
0ed6f095 761//___________________________________________________________________
b37d601d 762Bool_t AliTRDv0Info::Armenteros(const AliESDv0 *esdv0, Int_t decay){
0ed6f095 763 //
764 // computes the Armenteros variables for given V0
765 //
766 Double_t mn[3] = {0,0,0};
767 Double_t mp[3] = {0,0,0};
768 Double_t mm[3] = {0,0,0};
769
770 if(V0SignCheck()){
771 esdv0->GetNPxPyPz(mn[0],mn[1],mn[2]); //reconstructed cartesian momentum components of negative daughter
772 esdv0->GetPPxPyPz(mp[0],mp[1],mp[2]); //reconstructed cartesian momentum components of positive daughter
773 }
774 else{
775 esdv0->GetPPxPyPz(mn[0],mn[1],mn[2]); //reconstructed cartesian momentum components of negative daughter
776 esdv0->GetNPxPyPz(mp[0],mp[1],mp[2]); //reconstructed cartesian momentum components of positive daughter
777 }
778 esdv0->GetPxPyPz(mm[0],mm[1],mm[2]); //reconstructed cartesian momentum components of mother
779
780 TVector3 vecN(mn[0],mn[1],mn[2]);
781 TVector3 vecP(mp[0],mp[1],mp[2]);
782 TVector3 vecM(mm[0],mm[1],mm[2]);
783
784 Double_t thetaP = acos((vecP * vecM)/(vecP.Mag() * vecM.Mag()));
785 Double_t thetaN = acos((vecN * vecM)/(vecN.Mag() * vecM.Mag()));
786
787 Double_t alfa = ((vecP.Mag())*cos(thetaP)-(vecN.Mag())*cos(thetaN))/
788 ((vecP.Mag())*cos(thetaP)+(vecN.Mag())*cos(thetaN)) ;
789 Double_t qt = vecP.Mag()*sin(thetaP);
790
791 Float_t ap[2];
792 ap[0] = alfa;
793 ap[1] = qt;
794
64d57299 795 Double_t lCutAP[2];//Lambda/Anti-Lambda cuts
0ed6f095 796 if(decay == 0){
797 // armenteros cuts
798 const Double_t cutAlpha[2] = {0.35, 0.45}; // [0.35, 0.45]
799 const Double_t cutQT = 0.015;
800 if(TMath::Abs(ap[0]) > cutAlpha[0] && TMath::Abs(ap[0]) < cutAlpha[1]) return kFALSE;
801
802 if(ap[1] > cutQT) return kFALSE;
803 }
804
805 else if(decay == 1){
806 const Double_t cutQT = 0.1075;
807 const Double_t cutAP = 0.22 * TMath::Sqrt( TMath::Abs( (1-ap[0]*ap[0]/(0.92*0.92)) ) );
808 if(ap[1] < cutQT) return kFALSE;
809 if(ap[1] > cutAP) return kFALSE;
810 }
811 else if(decay == 2){
812 const Double_t cutQT = 0.03;
813 const Double_t cutAlpha = 0.7; // VERY strong - should supress the overlap with K0
64d57299 814 lCutAP[0] = 1.0 - (ap[0]-0.7 * ap[0]-0.7)*1.1 - 0.87;
0ed6f095 815 if(TMath::Abs(ap[0]) > cutAlpha) return kFALSE;
816 if(ap[1] < cutQT) return kFALSE;
64d57299 817 if(ap[1] > lCutAP[0]) return kFALSE;
0ed6f095 818
819 }
820 else if(decay == 3){
821 const Double_t cutQT = 0.03;
822 const Double_t cutAlpha = 0.7; // VERY strong - should supress the overlap with K0
64d57299 823 lCutAP[1] = 1.0 - (ap[0]+0.7 * ap[0]+0.7)*1.1 - 0.87;
0ed6f095 824 if(TMath::Abs(ap[0]) > cutAlpha) return kFALSE;
825 if(ap[1] < cutQT) return kFALSE;
64d57299 826 if(ap[1] > lCutAP[1]) return kFALSE;
0ed6f095 827 }
828 return kTRUE;
829}
5047978d 830
1ee39b3a 831//_________________________________________________
3d19c1b0 832Int_t AliTRDv0Info::GetPID(Int_t ipart, AliTRDtrackInfo *track)
833{
0ed6f095 834 // Decides if track is accepted for one of the reference data samples
5047978d 835
0ed6f095 836 Int_t cutCode = -99;
3d19c1b0 837 if(!(track)) {
838 AliError("No track info");
b9ddd472 839 return -1;
3d19c1b0 840 }
841 if(!HasTrack(track)){
842 AliDebug(2, "Track not attached to v0.");
0ed6f095 843 return -2;
3d19c1b0 844 }
3d19c1b0 845
846 //translate ipart to decay (Anti-Lambda will be treated separately)
1ee39b3a 847 Int_t iDecay = -1;
3d19c1b0 848 switch(ipart){
849 case AliPID::kElectron: iDecay = kGamma; break;
850 case AliPID::kPion: iDecay = kK0s; break;
851 case AliPID::kProton: iDecay = kLambda; break;
852 default:
b2f4ab8d 853 AliDebug(1, Form("Hypothesis \"ipart=%d\" not handled", ipart));
0ed6f095 854 return -3;
3d19c1b0 855 }
1ee39b3a 856
3d19c1b0 857 //... it fulfills our quality criteria
0ed6f095 858 if(!(fQuality == 1)) return -4;
3d19c1b0 859 //... distance of closest approach between daughters is reasonably small
0ed6f095 860 if((fDCA > fUpDCA[iDecay])) return -5;
3d19c1b0 861 //... pointing angle between momentum of mother particle and vector from prim. to sec. vertex is small
0ed6f095 862 if((fPointingAngle > fUpPointingAngle[iDecay])) return -6;
3d19c1b0 863 //... x-y plane distance of decay point to prim. vertex is bigger than a certain minimum value (for conversions)
0ed6f095 864 if((fRadius < fDownRadius[iDecay])) return -7;
3d19c1b0 865 //...or smaller than a maximum value (for K0s)
0ed6f095 866 if((fRadius > fUpRadius[iDecay])) return -8;
3d19c1b0 867 //... opening angle is close enough to zero (for conversions)
0ed6f095 868 if((fOpenAngle > fUpOpenAngle[iDecay])) return -9;
3d19c1b0 869 //... Psi-pair angle is close enough to zero(for conversions)
0ed6f095 870 if((TMath::Abs(fPsiPair) > fUpPsiPair[iDecay])) return -10;
871
3d19c1b0 872
873
874 //Mother momentum slots above/below 2.5 GeV
875 Int_t iPSlot(fV0Momentum > 2.5);
b9ddd472 876 Int_t trackID(track->GetTrackId());
3d19c1b0 877
878 //specific cut criteria :
879 if(ipart == AliPID::kProton) {
0ed6f095 880 if((fInvMass[kK0s] < fUpInvMass[kK0s][iPSlot]) && (fInvMass[kK0s] > fDownInvMass[kK0s])) return -11;//explicit exclusion of K0s decays
5047978d 881 if(fOpenAngle < (0.3 - 0.2*fV0Momentum)) return -9;
0ed6f095 882
3d19c1b0 883 //for proton sample: separate treatment of Lamba and Anti-Lambda decays:
884 //for Anti-Lambda:
885 //Combined PID likelihoods high enough for pi+ and anti-proton ; invariant mass calculated postulating these two particle species...
0ed6f095 886 //if((fComPID[kNeg][AliPID::kProton] > fDownComPIDneg[AliPID::kProton]) && (fComPID[kPos][AliPID::kPion] > fDownComPIDposPart[AliPID::kPion])) {
887 //if((fDetPID[kNeg][kTPC][AliPID::kProton] > fDownTPCPIDneg[AliPID::kProton]) && (fDetPID[kPos][kTPC][AliPID::kPion] > fDownTPCPIDpos[AliPID::kPion])){
888 if((TPCdEdxCuts(ipart, track))){//momentary solution: direct cut on TPC dE/dx
889 if(fNindex == trackID) {//we're only interested in the anti-proton
5047978d 890 if(fArmenteros[kAntiLambda]){//Armenteros condition has to be fulfilled
891 if(fChi2ndf[kAntiLambda] < fUpChi2ndf[kAntiLambda]){//Kalman filter Chi2/NDF not allowed to be too large
892 if((fInvMass[kAntiLambda] < fUpInvMass[kAntiLambda][iPSlot]) && (fInvMass[kAntiLambda] > fDownInvMass[kAntiLambda])){
893 return 1;
894 } else cutCode = -15;
895 } else cutCode =-14;
896 } else cutCode = -13;
3d19c1b0 897 }
5047978d 898 } else cutCode = -12;
3d19c1b0 899 //for Lambda:
900 //TPC PID likelihoods high enough for pi- and proton ; invariant mass calculated accordingly
0ed6f095 901 //if((fComPID[kNeg][AliPID::kPion] > fDownComPIDnegPart[AliPID::kPion]) && (fComPID[kPos][AliPID::kProton] > fDownComPIDpos[AliPID::kProton])) {
902 //if((fDetPID[kNeg][kTPC][AliPID::kPion] > fDownTPCPIDneg[AliPID::kPion]) && (fDetPID[kPos][kTPC][AliPID::kProton] > fDownTPCPIDpos[AliPID::kProton])){
903 if((TPCdEdxCuts(ipart, track))){//momentary solution: direct TPC dE/dx cuts
3d19c1b0 904 if(fPindex == trackID) {
5047978d 905 if(fArmenteros[kLambda]){
906 if(fChi2ndf[kLambda] < fUpChi2ndf[kLambda]){
907 if((fInvMass[kLambda] < fUpInvMass[kLambda][iPSlot]) && (fInvMass[kLambda] > fDownInvMass[kLambda])){
908 return 1;
909 } else cutCode = -15;
910 } else cutCode = -14;
911 } else cutCode = -13;
3d19c1b0 912 }
5047978d 913 } else cutCode = -12;
0ed6f095 914 return cutCode;
3d19c1b0 915 }
0ed6f095 916
3d19c1b0 917 //for K0s decays: equal TPC PID likelihood criteria for both daughters ; invariant mass calculated postulating two pions
918 if(ipart == AliPID::kPion) {
5047978d 919 if(fOpenAngle < (1.0/(fV0Momentum + 0.3) - 0.1)) return -9;
3d19c1b0 920 //explicit exclusion of Lambda decays
0ed6f095 921 if((fInvMass[kLambda] < fUpInvMass[kLambda][iPSlot]) && (fInvMass[kLambda] > fDownInvMass[kLambda])) return -11;
3d19c1b0 922 //explicit exclusion of Anti-Lambda decays
0ed6f095 923 if((fInvMass[kAntiLambda] < fUpInvMass[kAntiLambda][iPSlot]) && (fInvMass[kAntiLambda] > fDownInvMass[kAntiLambda])) return -11;
0ed6f095 924 //if((fDetPID[kNeg][kTPC][ipart] < fDownTPCPIDneg[ipart]) || (fDetPID[kPos][kTPC][ipart] < fDownTPCPIDpos[ipart])) return -12;
925 if(!(TPCdEdxCuts(ipart, track))){//momentary solution: direct TPC dE/dx cuts
926 return -12;
927 }
3d19c1b0 928 }
0ed6f095 929
930
3d19c1b0 931 //for photon conversions: equal combined PID likelihood criteria for both daughters ; invariant mass calculated postulating two electrons
932 //No Lambda/K0s exclusion is provided, since these contributions hardly ever interfere with gamma invariant mass!
0ed6f095 933 //Float_t momentum(track->GetESDinfo()->GetOuterParam()->P());
3d19c1b0 934 if(ipart == AliPID::kElectron) {
0ed6f095 935 //if(momentum > 1.75) {//since combined PID performs a little worse in simulations than TPC standalone for higher momenta, ONLY TPC PID is used here
936 //if((fDetPID[kNeg][kTPC][ipart] < fDownTPCPIDneg[ipart]) || (fDetPID[kPos][kTPC][ipart] < fDownTPCPIDpos[ipart])) return -12;
937 //} else {//for low momenta, combined PID from TOF and TPC is used to get rid of proton contamination
938 //if((fComPID[kNeg][ipart] > fDownComPIDneg[ipart]) && (fComPID[kPos][ipart] > fDownComPIDpos[ipart])) return 1;
939 //}
940 if(!(TPCdEdxCuts(ipart, track))){//momentary solution for direct TPC dE/dx cut
941 return -12;
1ee39b3a 942 }
0ed6f095 943 }
944
0ed6f095 945 //Armenteros-Polanski cut
946 if(!(fArmenteros[iDecay])) return -13;
0ed6f095 947 //Kalman filter Chi2/NDF cut
948 if(fChi2ndf[iDecay] > fUpChi2ndf[iDecay]) return -14;
0ed6f095 949 //Invariant mass cut for K0s and photons, assuming two pions/two electrons as daughters:
5047978d 950 if((fInvMass[iDecay] > fUpInvMass[iDecay][iPSlot]) || (fInvMass[iDecay] < fDownInvMass[iDecay])) return -15;
0ed6f095 951
952 return 1;
1ee39b3a 953}
3d19c1b0 954
955
1ee39b3a 956//_________________________________________________
b9ddd472 957void AliTRDv0Info::Print(Option_t *opt) const
1ee39b3a 958{
64d57299 959 //prints text for debugging etc.
5047978d 960 printf("V0 legs :: %4d[+] %4d[-]\n", fPindex, fNindex);
961 printf(" Decay :: Gamma[%c] K0s[%c] Lambda[%c] AntiLambda[%c]\n",
962 IsDecay(kGamma)?'y':'n', IsDecay(kK0s)?'y':'n', IsDecay(kLambda)?'y':'n', IsDecay(kAntiLambda)?'y':'n');
963 printf(" Kine :: DCA[%5.3f] Radius[%5.3f]\n", fDCA, fRadius);
964 printf(" Angle :: Pointing[%5.3f] Open[%5.3f] Psi[%5.3f]\n", fPointingAngle, fOpenAngle, fPsiPair);
b9ddd472 965 if(strcmp(opt, "a")!=0) return;
5047978d 966 printf(" PID ::\n"
967 " ITS TPC TOF COM\n");
b9ddd472 968 for(Int_t idt=0; idt<kNDaughters; idt++){
b9ddd472 969 for(Int_t is(0); is<AliPID::kSPECIES; is++){
5047978d 970 printf(" %s%c%s", AliPID::ParticleShortName(is), idt?'-':'+', (is==1||is==2)?" ":" ");
b9ddd472 971 for(Int_t id(0); id<kNDetectors; id++){
972 printf("%5.1f ", 1.e2*fDetPID[idt][id][is]);
973 }
974 printf("%5.1f\n", 1.e2*fComPID[idt][is]);
975 }
976 }
977}
1ee39b3a 978
b9ddd472 979//_________________________________________________
980void AliTRDv0Info::SetV0tracks(AliESDtrack *p, AliESDtrack *n)
981{
64d57299 982 //sets the two daughter trex and their indices
b9ddd472 983 fTrackP = p; fPindex = p->GetID();
984 fTrackN = n; fNindex = n->GetID();
1ee39b3a 985}
64d57299 986//_________________________________________________
b9ddd472 987
64d57299 988AliESDtrack *AliTRDv0Info::GetV0Daughter(Int_t sign)
989{
990 //Gets positive of negative daughter of decay
991 if(sign>0)
992 return fTrackP;
993 else if(sign < 0)
994 return fTrackN;
b9ddd472 995
64d57299 996 return 0;
997}