]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG3/hfe/AliHFEV0cuts.cxx
Excluding lowest pt point for K and p, ratios fit/data are now shown in 3 canvases
[u/mrichter/AliRoot.git] / PWG3 / hfe / AliHFEV0cuts.cxx
CommitLineData
c04c80e6 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// * 20/04/2010 *
16// Class for optimising and applying V0 cuts to obtain clean V0 samples
17// Compatible with ESDs only
18//
19// Authors:
20// Matus Kalisky <matus.kalisky@cern.ch>
21//
22
23#include "TDatabasePDG.h"
24
25#include "AliESDtrack.h"
26#include "AliMCEvent.h"
27#include "AliESDv0.h"
28#include "AliKFParticle.h"
29#include "AliKFVertex.h"
30
31#include "AliHFEcollection.h"
32
33#include "AliHFEV0cuts.h"
34
35ClassImp(AliHFEV0cuts)
36
37//________________________________________________________________
38AliHFEV0cuts::AliHFEV0cuts():
39 fQA(NULL)
40 , fMCEvent(NULL)
41 , fInputEvent(NULL)
42 , fPrimaryVertex(NULL)
43 , fIsMC(kFALSE)
44{
45
46 //
47 // Default constructor
48 //
49
50
51}
52//________________________________________________________________
53AliHFEV0cuts::~AliHFEV0cuts()
54{
55 //
56 // destructor
57 //
58 if (fQA) delete fQA;
59}
60
61//________________________________________________________________
62AliHFEV0cuts::AliHFEV0cuts(const AliHFEV0cuts &ref):
63 TObject(ref)
64 , fQA(NULL)
65 , fMCEvent(NULL)
66 , fInputEvent(NULL)
67 , fPrimaryVertex(NULL)
68 , fIsMC(kFALSE)
69{
70 //
71 // Copy constructor
72 //
73 ref.Copy(*this);
74}
75//________________________________________________________________
76AliHFEV0cuts &AliHFEV0cuts::operator=(const AliHFEV0cuts &ref){
77 //
78 // Assignment operator
79 //
80 if(this != &ref)
81 ref.Copy(*this);
82 return *this;
83}
84//________________________________________________________________
85void AliHFEV0cuts::Copy(TObject &ref) const{
86 //
87 // Copy function
88 //
89 AliHFEV0cuts &target = dynamic_cast<AliHFEV0cuts &>(ref);
90
91 if(fQA) target.fQA = dynamic_cast<AliHFEcollection *>(fQA->Clone());
92
93 if(target.fMCEvent) delete target.fMCEvent;
94 target.fMCEvent = new AliMCEvent;
95
96 if(target.fPrimaryVertex) delete target.fPrimaryVertex;
97 target.fPrimaryVertex = new AliKFVertex;
98
99 TObject::Copy(ref);
100
101}
102//___________________________________________________________________
103void AliHFEV0cuts::Init(const char* name){
104 //
105 // initialize the output objects and create histograms
106 //
107
108 //
109 // all the "h_cut_XXX" histograms hare cut value distributions:
110 // [0] for all candidates
111 // [1] jus before the cut on given variable was applied, but after all the previous cuts
112 //
113
114 fQA = new AliHFEcollection("fQA", name);
115
116
117 // common for all V0s
118 fQA->CreateTH2Fvector1(2, "h_all_AP", "armenteros plot for all V0 candidates", 200, -1, 1, 200, 0, 0.25);
119
120 // gammas
121 fQA->CreateTH1Fvector1(2, "h_cut_Gamma_CosPoint", "Gamma Cosine pointing angle; cos point. angle; counts", 100, 0, 0.1);
122 fQA->CreateTH1Fvector1(2, "h_cut_Gamma_DCA", "DCA between the gamma daughters; dca (cm); counts", 100, 0, 2);
123 fQA->CreateTH1Fvector1(2, "h_cut_Gamma_VtxR", "Radius of the gamma conversion vertex; r (cm); counts", 1000, 0, 100);
124 fQA->CreateTH1Fvector1(2, "h_cut_Gamma_OA", "opening angle of the gamma products; opening angle (rad); counts", 100, 0, 1);
125 fQA->CreateTH1Fvector1(2, "h_cut_Gamma_PP", "gamma psi pair angle; psi pairangle (rad); counts", 100, 0, 2);
126 fQA->CreateTH1Fvector1(2, "h_cut_Gamma_Chi2", "gamma Chi2/NDF; Chi2/NDF; counts", 100, 0, 25);
127 fQA->CreateTH1Fvector1(9, "h_Gamma_Mass", "Invariant mass of gammas; mass (GeV/c^{2}); counts", 100, 0, 0.2);
128
129
130 // kaons
131 fQA->CreateTH1Fvector1(2, "h_cut_K0_CosPoint", "K0 Cosine pointing angle; cos point. angle; counts", 100, 0, 0.1);
132 fQA->CreateTH1Fvector1(2, "h_cut_K0_DCA", "DCA between the K0 daughters; dca (cm); counts", 100, 0, 2);
133 fQA->CreateTH1Fvector1(2, "h_cut_K0_VtxR", "Radius of the K0 decay vertex; r (cm); counts", 1000, 0, 100);
134 fQA->CreateTH1Fvector1(2, "h_cut_K0_Chi2", "K0 Chi2/NDF; Chi2/NDF; counts", 100, 0, 25);
135 fQA->CreateTH2Fvector1(2, "h_cut_K0_AP", "Armenteros plot for K0 candidates; #alpha; q_{T} (GeV/c)", 100, -1, 1, 100, 0, 0.3);
136 fQA->CreateTH1Fvector1(8, "h_K0_Mass", "Invariant mass of K0; mass (GeV/c^{2}); counts", 125, 0.45, 0.55);
137
138 // lambda
139 fQA->CreateTH1Fvector1(2, "h_cut_L_CosPoint", "L Cosine pointing angle; cos point. angle; counts", 100, 0, 0.1);
140 fQA->CreateTH1Fvector1(2, "h_cut_L_DCA", "DCA between the L daughters; dca (cm); counts", 100, 0, 2);
141 fQA->CreateTH1Fvector1(2, "h_cut_L_VtxR", "Radius of the L decay vertex; r (cm); counts", 1000, 0, 100);
142 fQA->CreateTH1Fvector1(2, "h_cut_L_Chi2", "L Chi2/NDF; Chi2/NDF; counts", 100, 0, 25);
143 fQA->CreateTH2Fvector1(2, "h_cut_L_AP", "Armenteros plot for Lambda candidates; #alpha; q_{T} (GeV/c)", 100, -1, 1, 100, 0, 0.3);
144 fQA->CreateTH2Fvector1(2, "h_cut_AL_AP", "Armenteros plot for anti Lambda candidates; #alpha; q_{T} (GeV/c)", 100, -1, 1, 100, 0, 0.3);
145 fQA->CreateTH2Fvector1(2, "h_cut_Gamma_AP", "Armenteros plot for gamma candidates; #alpha; q_{T} (GeV/c)", 100, -1, 1, 100, 0, 0.3);
146 fQA->CreateTH1Fvector1(9, "h_L_Mass", "Invariant mass of L; mass (GeV/c^{2}); counts", 60, 1.1, 1.13);
147 fQA->CreateTH1Fvector1(9, "h_AL_Mass", "Invariant mass of anti L; mass (GeV/c^{2}); counts", 60, 1.1, 1.13);
148
149 fQA->CreateTH2F("h_L_checks", "Lambda candidate check[0] -v- check[1]; check[0]; check[1]", 5, -0.75, 1.75, 6, -0.75, 1.75 );
150
151 // electrons
152 fQA->CreateTH1Fvector1(9, "h_Electron_P", "Momenta of conversion electrons -cuts-; P (GeV/c); counts", 50, 0.1, 20, 0);
153
154 // K0 pions
155 fQA->CreateTH1Fvector1(8, "h_PionK0_P", "Momenta of K0 pions -cuts-; P (GeV/c) counts;", 50, 0.1, 20, 0);
156
157 // L pions
158 fQA->CreateTH1Fvector1(9, "h_PionL_P", "Momenta of L pions -cuts-; P (GeV/c) counts;", 50, 0.1, 50, 0);
159
160 // L protons
161 fQA->CreateTH1Fvector1(9, "h_ProtonL_P", "Momenta of L protons -cuts-; P (GeV/c) counts;", 50, 0.1, 20, 0);
162
163 // single track cuts
164 fQA->CreateTH1F("h_ST_NclsTPC", "Number of TPC clusters", 161, -1, 160);
165 fQA->CreateTH1F("h_ST_TPCrefit", "TPC refit", 2, -0.5, 1.5);
166 fQA->CreateTH1F("h_ST_chi2TPCcls", "chi2 per TPC cluster", 100, 0, 10);
167 fQA->CreateTH1F("h_ST_TPCclsR", "TPC cluster ratio", 120, -0.1, 1.1);
168 fQA->CreateTH1F("h_ST_kinks", "kinks", 2, -0.5, 1.5);
169 fQA->CreateTH1F("h_ST_pt", "track pt", 100, 0.1, 20, 0);
170 fQA->CreateTH1F("h_ST_eta", "track eta", 100, -1.5, 1.5);
171
172 //
173 // possibly new cuts
174 //
175
176 // Gamma
177 fQA->CreateTH2Fvector1(2, "h_cut_Gamma_OAvP", "open. ang. of the Gamma daughters versus Gamma mom; Gamma p (GeV/c); opening angle (pions) (rad)", 100, 0.1, 10, 200, 0., 0.2);
178 // K0
179 fQA->CreateTH2Fvector1(2, "h_cut_K0_OAvP", "open. ang. of the K0 daughters versus K0 momentum; K0 p (GeV/c); opening angle (pions) (rad)", 100, 0.1, 10, 100, 0, 3.5);
180 // Lambda
181 fQA->CreateTH2Fvector1(2, "h_cut_L_OAvP", "open. ang. of the L daughters versus L momentum; Lambda p (GeV/c); openeing angle pion-proton (rad)", 100, 0.1, 10, 100, 0, 3.5);
182 fQA->CreateTH2Fvector1(2, "h_cut_L_rdp_v_mp", "relative L daughter mom -v- mother mom; L mom (GeV/c); relative daughter mom p2/p1", 100, 0.1, 10, 100, 0, 1);
183 fQA->CreateTH2Fvector1(2, "h_cut_L_qT_v_mp", "A-P q_{T} -v- Lambda momentum; mom (GeV/c); q_{T} GeV/c", 100, 0.1, 10, 50, 0., 0.12);
184
185
186 // THnSparse histograms
187
188 // THnSparse for the K0 mass
189 // to be looked at after merging run by run
190 // axes: mass, pt, theta, phi
191 {
192 Int_t nBin[4] = {100, 10, 10, 18};
193 Double_t nMin[4] = {0.45, 0.1, 0., 0.};
194 Double_t nMax[4] = {0.55, 10., TMath::Pi(), 2*TMath::Pi()};
195 TString htitle = "K0 sparse; mass (GeV/c^{2}); p_{T} (GeV/c); theta (rad); phi(rad)";
196 fQA->CreateTHnSparse("hK0", htitle, 4, nBin, nMin, nMax);
197 fQA->BinLogAxis("hK0", 1);
198 }
199
200}
201//________________________________________________________________
202Bool_t AliHFEV0cuts::TrackCutsCommon(AliESDtrack* track){
203 //
204 // singe track cuts commom for all particle candidates
205 //
206
207 if(!track) return kFALSE;
208
209 // status word
210 ULong_t status = track->GetStatus();
211
212
213 // No. of TPC clusters
214 fQA->Fill("h_ST_NclsTPC", track->GetTPCNcls());
215 if(track->GetTPCNcls() < 80) return kFALSE;
216
217 // TPC refit
218 if((status & AliESDtrack::kTPCrefit)){
219 fQA->Fill("h_ST_TPCrefit", 1);
220 }
221 if(!(status & AliESDtrack::kTPCrefit)){
222 fQA->Fill("h_ST_TPCrefit", 0);
223 return kFALSE;
224 }
225
226 // Chi2 per TPC cluster
227 Int_t nTPCclusters = track->GetTPCclusters(0);
228 Float_t chi2perTPCcluster = track->GetTPCchi2()/Float_t(nTPCclusters);
229 fQA->Fill("h_ST_chi2TPCcls", chi2perTPCcluster);
230 if(chi2perTPCcluster > 3.5) return kFALSE;
231
232 // TPC cluster ratio
233 Float_t cRatioTPC = track->GetTPCNclsF() > 0. ? static_cast<Float_t>(track->GetTPCNcls())/static_cast<Float_t> (track->GetTPCNclsF()) : 1.;
234 fQA->Fill("h_ST_TPCclsR", cRatioTPC);
235 if(cRatioTPC < 0.6) return kFALSE;
236
237 // kinks
238 fQA->Fill("h_ST_kinks", track->GetKinkIndex(0));
239 if(track->GetKinkIndex(0) != 0) return kFALSE;
240
241 // pt
242 fQA->Fill("h_ST_pt",track->Pt());
243 if(track->Pt() < 0.1 || track->Pt() > 100) return kFALSE;
244
245 // eta
246 fQA->Fill("h_ST_eta", track->Eta());
247 //if(TMath::Abs(track->Eta()) > 0.9) return kFALSE;
248
249 return kTRUE;
250}
251//________________________________________________________________
252Bool_t AliHFEV0cuts::V0CutsCommon(AliESDv0 *v0){
253 //
254 // V0 cuts common to all V0s
255 //
256
257 AliESDtrack* dN, *dP;
258
259 dP = dynamic_cast<AliESDtrack *>(fInputEvent->GetTrack(v0->GetPindex()));
260 dN = dynamic_cast<AliESDtrack *>(fInputEvent->GetTrack(v0->GetNindex()));
261
262 if(!dN || !dP) return kFALSE;
263
264 Int_t qP = dP->Charge();
265 Int_t qN = dN->Charge();
266
267 if((qP*qN) != -1) return kFALSE;
268
269 return kTRUE;
270}
271//________________________________________________________________
272Bool_t AliHFEV0cuts::GammaCuts(AliESDv0 *v0){
273 //
274 // gamma cuts
275 //
276
277 if(!v0) return kFALSE;
278
279 // loose cuts first
280 if(LooseRejectK0(v0) || LooseRejectLambda(v0)) return kFALSE;
281
282 AliVTrack* daughter[2];
283 Int_t pIndex = 0, nIndex = 0;
284 if(CheckSigns(v0)){
285 pIndex = v0->GetPindex();
286 nIndex = v0->GetNindex();
287 }
288 else{
289 pIndex = v0->GetNindex();
290 nIndex = v0->GetPindex();
291 }
292 daughter[0] = dynamic_cast<AliVTrack *>(fInputEvent->GetTrack(pIndex));
293 daughter[1] = dynamic_cast<AliVTrack *>(fInputEvent->GetTrack(nIndex));
294 if(!daughter[0] || !daughter[1]) return kFALSE;
295
296 AliKFParticle *kfMother = CreateMotherParticle(daughter[0], daughter[1], TMath::Abs(kElectron), TMath::Abs(kElectron));
297 if(!kfMother) return kFALSE;
298
299 // production vertex is set in the 'CreateMotherParticle' function
300 kfMother->SetMassConstraint(0, 0.001);
301
302 AliESDtrack* d[2];
303 d[0] = dynamic_cast<AliESDtrack*>(fInputEvent->GetTrack(pIndex));
304 d[1] = dynamic_cast<AliESDtrack*>(fInputEvent->GetTrack(nIndex));
305
306 Float_t iMass = v0->GetEffMass(0, 0);
307 Float_t iP = v0->P();
308 Float_t p[2] = {d[0]->GetP(), d[1]->GetP()};
309
310 // Armenteros
311 Float_t ap[2];
312 Armenteros(v0, ap);
313
314 // Cut values
315 const Double_t cutCosPoint[2] = {0., 0.03}; // ORG [0., 0.03]
316 const Double_t cutDCA[2] = {0., 0.25}; // ORG [0., 0.25]
317 const Double_t cutProdVtxR[2] = {6., 50.}; // ORG [6., 9999]
318 const Double_t cutOAngle[2] = {0, 0.1}; // ORG [0., 0.1]
319 const Double_t cutPsiPair[2] = {0., 0.05}; // ORG [0. 0.05]
320 const Double_t cutMass = 0.05; // ORG [0.05]
321 const Double_t cutChi2NDF = 7.; // ORG [7.]
322 // armenteros cuts
323 const Double_t cutAlpha[2] = {0.35, 0.45}; // [0.35, 0.45]
324 const Double_t cutQT = 0.015;
325
326 // Values
327
328 // cos pointing angle
329 Double_t cosPoint = v0->GetV0CosineOfPointingAngle();
330 cosPoint = TMath::ACos(cosPoint);
331
332 // DCA between daughters
333 Double_t dca = v0->GetDcaV0Daughters();
334
335 // Production vertex
336 Double_t x, y, z;
337 v0->GetXYZ(x,y,z);
338 Double_t r = TMath::Sqrt(x*x + y*y);
339
340 // Opening angle
341 Double_t oAngle = OpenAngle(v0);
342
343 // psi pair
344 Double_t psiPair = PsiPair(v0);
345
346 // V0 chi2/ndf
347 Double_t chi2ndf = kfMother->GetChi2()/kfMother->GetNDF();
348
349 if(kfMother) delete kfMother;
350
351 //
352 // Apply the cuts, produce QA plots (with mass cut)
353 //
354 fQA->Fill("h_Gamma_Mass", 0, iMass);
355 if(iMass < cutMass){
356 fQA->Fill("h_all_AP", 0, ap[0], ap[1]);
357 fQA->Fill("h_Electron_P", 0, p[0]);
358 fQA->Fill("h_Electron_P", 0, p[1]);
359 fQA->Fill("h_cut_Gamma_CosPoint", 0, cosPoint);
360 fQA->Fill("h_cut_Gamma_CosPoint", 1, cosPoint);
361 fQA->Fill("h_cut_Gamma_DCA", 0, dca);
362 fQA->Fill("h_cut_Gamma_VtxR", 0, r);
363 fQA->Fill("h_cut_Gamma_OA", 0, oAngle);
364 fQA->Fill("h_cut_Gamma_PP", 0, psiPair);
365 fQA->Fill("h_cut_Gamma_Chi2", 0, chi2ndf);
366 fQA->Fill("h_cut_Gamma_OAvP", 0, iP, oAngle);
367 fQA->Fill("h_cut_Gamma_AP", 0, ap[0], ap[1]);
368 }
369
370 if(cosPoint < cutCosPoint[0] || cosPoint > cutCosPoint[1]) return kFALSE;
371 fQA->Fill("h_Gamma_Mass", 1, iMass);
372 if(iMass < cutMass){
373 fQA->Fill("h_Electron_P", 1, p[0]);
374 fQA->Fill("h_Electron_P", 1, p[1]);
375 fQA->Fill("h_cut_Gamma_DCA", 1, dca);
376 }
377
378 if(dca < cutDCA[0] || dca > cutDCA[1]) return kFALSE;
379 fQA->Fill("h_Gamma_Mass", 2, iMass);
380 if(iMass < cutMass){
381 fQA->Fill("h_Electron_P", 2, p[0]);
382 fQA->Fill("h_Electron_P", 2, p[1]);
383 fQA->Fill("h_cut_Gamma_VtxR", 1, r);
384 }
385
386 if(r < cutProdVtxR[0] || r > cutProdVtxR[1]) return kFALSE;
387 fQA->Fill("h_Gamma_Mass", 3, iMass);
388 if(iMass < cutMass){
389 fQA->Fill("h_Electron_P", 3, p[0]);
390 fQA->Fill("h_Electron_P", 3, p[1]);
391 fQA->Fill("h_cut_Gamma_OA", 1, oAngle);
392 }
393
394 if(oAngle < cutOAngle[0] || oAngle > cutOAngle[1]) return kFALSE;
395 fQA->Fill("h_Gamma_Mass", 4, iMass);
396 if(iMass < cutMass){
397 fQA->Fill("h_Electron_P", 4, p[0]);
398 fQA->Fill("h_Electron_P", 4, p[1]);
399 fQA->Fill("h_cut_Gamma_PP", 1, psiPair);
400 }
401
402 if(psiPair < cutPsiPair[0] || psiPair > cutPsiPair[1]) return kFALSE;
403 fQA->Fill("h_Gamma_Mass", 5, iMass);
404 if(iMass < cutMass){
405 fQA->Fill("h_Electron_P", 5, p[0]);
406 fQA->Fill("h_Electron_P", 5, p[1]);
407 fQA->Fill("h_cut_Gamma_Chi2", 1, chi2ndf);
408 }
409
410 if(chi2ndf > cutChi2NDF) return kFALSE;
411 fQA->Fill("h_Gamma_Mass", 6, iMass);
412 if(iMass < cutMass){
413 fQA->Fill("h_Electron_P", 6, p[0]);
414 fQA->Fill("h_Electron_P", 6, p[1]);
415 fQA->Fill("h_cut_Gamma_OAvP", 1, iP, oAngle);
416 fQA->Fill("h_all_AP", 1, ap[0], ap[1]);
417 fQA->Fill("h_cut_Gamma_AP", 1, ap[0], ap[1]);
418 }
419
420 if(TMath::Abs(ap[0]) > cutAlpha[0] && TMath::Abs(ap[0]) < cutAlpha[1]) return kFALSE;
421 fQA->Fill("h_Gamma_Mass", 7, iMass);
422 if(iMass < cutMass){
423 fQA->Fill("h_Electron_P", 7, p[0]);
424 fQA->Fill("h_Electron_P", 7, p[1]);
425 }
426
427 if(ap[1] > cutQT) return kFALSE;
428 fQA->Fill("h_Gamma_Mass", 8, iMass);
429 if(iMass < cutMass){
430 fQA->Fill("h_Electron_P", 8, p[0]);
431 fQA->Fill("h_Electron_P", 8, p[1]);
432 }
433
434
435 if(iMass > cutMass) return kFALSE;
436
437 // all cuts passed
438
439 return kTRUE;
440}
441//________________________________________________________________
442Bool_t AliHFEV0cuts::K0Cuts(AliESDv0 *v0){
443 //
444 // K0 cuts
445 //
446
447 if(!v0) return kFALSE;
448
449 const Double_t cK0mass=TDatabasePDG::Instance()->GetParticle(kK0Short)->Mass(); // PDG K0s mass
450 AliVTrack* daughter[2];
451 Int_t pIndex = 0, nIndex = 0;
452 if(CheckSigns(v0)){
453 pIndex = v0->GetPindex();
454 nIndex = v0->GetNindex();
455 }
456 else{
457 pIndex = v0->GetNindex();
458 nIndex = v0->GetPindex();
459 }
460
461 daughter[0] = dynamic_cast<AliVTrack *>(fInputEvent->GetTrack(pIndex));
462 daughter[1] = dynamic_cast<AliVTrack *>(fInputEvent->GetTrack(nIndex));
463 if(!daughter[0] || !daughter[1]) return kFALSE;
464
465 AliKFParticle *kfMother = CreateMotherParticle(daughter[0], daughter[1], TMath::Abs(kPiPlus), TMath::Abs(kPiPlus));
466 if(!kfMother) return kFALSE;
467 // production vertex is set in the 'CreateMotherParticle' function
468 kfMother->SetMassConstraint(cK0mass, 0.);
469
470 AliESDtrack* d[2];
471 d[0] = dynamic_cast<AliESDtrack*>(fInputEvent->GetTrack(pIndex));
472 d[1] = dynamic_cast<AliESDtrack*>(fInputEvent->GetTrack(nIndex));
473
474 Float_t iMass = v0->GetEffMass(2, 2);
475 Float_t iP = v0->P();
476 Float_t p[2] = {d[0]->GetP(), d[1]->GetP()};
477 Double_t theta = v0->Theta();
478 Double_t phi = v0->Phi();
479 Double_t pt = v0->Pt();
480 Double_t data[4] = {0., 0., 0., 0.};
481 // Cut values
482 const Double_t cutCosPoint[2] = {0., 0.03}; // ORG [0., 0.03]
483 const Double_t cutDCA[2] = {0., 0.2}; // ORG [0., 0.1]
484 const Double_t cutProdVtxR[2] = {2.0, 30.}; // ORG [0., 8.1]
485 const Double_t cutMass[2] = {0.49, 0.51}; // ORG [0.485, 0.51]
486 const Double_t cutChi2NDF = 5.; // ORG [7.]
487 const Double_t cutOAngleP = (1.0/(iP + 0.3) - 0.1); // momentum dependent min. OAngle ~ 1/x
488 // cundidate cuts
489 // armenteros plot
490 const Double_t cutQT = 0.1075;
491 // elipse cut - see bellow
492
493 // Values
494
495 // cos pointing angle
496 Double_t cosPoint = v0->GetV0CosineOfPointingAngle();
497 cosPoint = TMath::ACos(cosPoint);
498
499 // DCA between daughters
500 Double_t dca = v0->GetDcaV0Daughters();
501
502 // Production vertex
503 Double_t x, y, z;
504 v0->GetXYZ(x,y,z);
505
506 Double_t r = TMath::Sqrt(x*x + y*y);
507
508 // V0 chi2/ndf
509 Double_t chi2ndf = kfMother->GetChi2()/kfMother->GetNDF();
510
511 if(kfMother) delete kfMother;
512
513 // Opening angle
514 Double_t oAngle = OpenAngle(v0);
515
516 // Armenteros
517 Float_t ap[2];
518 Armenteros(v0, ap);
519 const Double_t cutAP = 0.22 * TMath::Sqrt( TMath::Abs( (1-ap[0]*ap[0]/(0.92*0.92)) ) );
520
521
522 //
523 // Apply the cuts, produce QA plots (with mass cut)
524 //
525
526 fQA->Fill("h_K0_Mass", 0, iMass);
527 if(iMass > cutMass[0] && iMass < cutMass[1]){
528 fQA->Fill("h_all_AP", 0, ap[0], ap[1]);
529 fQA->Fill("h_PionK0_P", 0, p[0]);
530 fQA->Fill("h_PionK0_P", 0, p[1]);
531 fQA->Fill("h_cut_K0_CosPoint", 0, cosPoint);
532 fQA->Fill("h_cut_K0_CosPoint", 1, cosPoint);
533 fQA->Fill("h_cut_K0_DCA", 0, dca);
534 fQA->Fill("h_cut_K0_VtxR", 0, r);
535 fQA->Fill("h_cut_K0_Chi2", 0, chi2ndf);
536 fQA->Fill("h_cut_K0_OAvP", 0, iP, oAngle);
537 fQA->Fill("h_cut_K0_AP", 0, ap[0], ap[1]);
538 }
539
540 if(cosPoint < cutCosPoint[0] || cosPoint > cutCosPoint[1]) return kFALSE;
541 fQA->Fill("h_K0_Mass", 1, iMass);
542 if(iMass > cutMass[0] && iMass < cutMass[1]){
543 fQA->Fill("h_PionK0_P", 1, p[0]);
544 fQA->Fill("h_PionK0_P", 1, p[1]);
545 fQA->Fill("h_cut_K0_DCA", 1, dca);
546 }
547
548 if(dca < cutDCA[0] || dca > cutDCA[1]) return kFALSE;
549 fQA->Fill("h_K0_Mass", 2, iMass);
550 if(iMass > cutMass[0] && iMass < cutMass[1]){
551 fQA->Fill("h_PionK0_P", 2, p[0]);
552 fQA->Fill("h_PionK0_P", 2, p[1]);
553 fQA->Fill("h_cut_K0_VtxR", 1, r);
554 }
555
556
557 if(r < cutProdVtxR[0] || r > cutProdVtxR[1]) return kFALSE;
558 fQA->Fill("h_K0_Mass", 3, iMass);
559 if(iMass > cutMass[0] && iMass < cutMass[1]){
560 fQA->Fill("h_PionK0_P", 3, p[0]);
561 fQA->Fill("h_PionK0_P", 3, p[1]);
562 fQA->Fill("h_cut_K0_Chi2", 1, chi2ndf);
563 }
564
565 if(chi2ndf > cutChi2NDF) return kFALSE;
566 fQA->Fill("h_K0_Mass", 4, iMass);
567 if(iMass > cutMass[0] && iMass < cutMass[1]){
568 fQA->Fill("h_PionK0_P", 4, p[0]);
569 fQA->Fill("h_PionK0_P", 4, p[1]);
570 fQA->Fill("h_cut_K0_OAvP", 1, iP, oAngle);
571 }
572
573 if(oAngle < cutOAngleP) return kFALSE;
574 fQA->Fill("h_K0_Mass", 5, iMass);
575 if(iMass > cutMass[0] && iMass < cutMass[1]){
576 fQA->Fill("h_PionK0_P", 5, p[0]);
577 fQA->Fill("h_PionK0_P", 5, p[1]);
578 fQA->Fill("h_cut_K0_AP", 1, ap[0], ap[1]);
579 fQA->Fill("h_all_AP", 1, ap[0], ap[1]);
580 }
581
582 if(ap[1] < cutQT) return kFALSE;
583 fQA->Fill("h_K0_Mass", 6, iMass);
584 if(iMass > cutMass[0] && iMass < cutMass[1]){
585 fQA->Fill("h_PionK0_P", 6, p[0]);
586 fQA->Fill("h_PionK0_P", 6, p[1]);
587 }
588
589 if(ap[1] > cutAP) return kFALSE;
590 fQA->Fill("h_K0_Mass", 7, iMass);
591 if(iMass > cutMass[0] && iMass < cutMass[1]){
592 fQA->Fill("h_PionK0_P", 7, p[0]);
593 fQA->Fill("h_PionK0_P", 7, p[1]);
594 }
595
596 data[0] = iMass;
597 data[1] = pt;
598 data[2] = theta;
599 data[3] = phi;
600 //printf("-D: m: %f, pT: %f, theta: %f, phi: %f\n", invMass, mPt, theta, phi);
601 fQA->Fill("hK0", data);
602
603
604 if(iMass < cutMass[0] || iMass > cutMass[1]) return kFALSE;
605
606 // all cuts passed
607
608 return kTRUE;
609}
610//________________________________________________________________
611Bool_t AliHFEV0cuts::LambdaCuts(AliESDv0 *v0, Bool_t &isLambda ){
612 //
613 // Lambda cuts - decision on Lambda - AntiLambda is taken too
614 //
615 // discrimination between lambda and antilambda - correlation of the following variables necessary:
616 // - momentum of the proton AND momentum of the pion (proton momentum is allways larger)
617 // - mass of the mother particle
618
619 if(!v0) return kFALSE;
620
621 // loose cuts first
622 if(LooseRejectK0(v0) || LooseRejectGamma(v0)) return kFALSE;
623
624 const Double_t cL0mass=TDatabasePDG::Instance()->GetParticle(kLambda0)->Mass(); // PDG lambda mass
625
626 AliVTrack* daughter[2];
627 Int_t pIndex = 0, nIndex = 0;
628 Float_t mMass[2] = {-1., -1.};
629 if(CheckSigns(v0)){
630 pIndex = v0->GetPindex();
631 nIndex = v0->GetNindex();
632 mMass[0] = v0->GetEffMass(4, 2);
633 mMass[1] = v0->GetEffMass(2, 4);
634 }
635 else{
636 pIndex = v0->GetNindex();
637 nIndex = v0->GetPindex();
638 mMass[0] = v0->GetEffMass(2, 4);
639 mMass[1] = v0->GetEffMass(4, 2);
640 }
641
642 daughter[0] = dynamic_cast<AliVTrack *>(fInputEvent->GetTrack(pIndex));
643 daughter[1] = dynamic_cast<AliVTrack *>(fInputEvent->GetTrack(nIndex));
644 if(!daughter[0] || !daughter[1]) return kFALSE;
645
646 AliKFParticle *kfMother[2] = {0x0, 0x0};
647 // Lambda
648 kfMother[0] = CreateMotherParticle(daughter[0], daughter[1], TMath::Abs(kProton), TMath::Abs(kPiPlus));
649 if(!kfMother[0]) return kFALSE;
650
651 // production vertex is set in the 'CreateMotherParticle' function
652 kfMother[0]->SetMassConstraint(cL0mass, 0.);
653
654 // Anti Lambda
655 kfMother[1] = CreateMotherParticle(daughter[0], daughter[1], TMath::Abs(kPiPlus), TMath::Abs(kProton));
656 if(!kfMother[1]) return kFALSE;
657 // production vertex is set in the 'CreateMotherParticle' function
658 kfMother[1]->SetMassConstraint(cL0mass, 0.);
659
660 Float_t dMass[2] = {TMath::Abs(mMass[0] - cL0mass), TMath::Abs(mMass[1] - cL0mass)};
661
662 AliESDtrack* d[2];
663 d[0] = dynamic_cast<AliESDtrack*>(fInputEvent->GetTrack(pIndex));
664 d[1] = dynamic_cast<AliESDtrack*>(fInputEvent->GetTrack(nIndex));
665 if(!d[0] || !d[1]) return kFALSE;
666
667 Float_t p[2] = {d[0]->GetP(), d[1]->GetP()};
668
669 // check the 3 lambda - antilambda variables
670 Int_t check[2] = {-1, -1}; // 0 : lambda, 1 : antilambda
671 // 1) momentum of the daughter particles - proton is expected to have higher momentum than pion
672 check[0] = (p[0] > p[1]) ? 0 : 1;
673 // 2) mass of the mother particle
674 check[1] = (dMass[0] < dMass[1]) ? 0 : 1;
675 fQA->Fill("h_L_checks", check[0]*1.0, check[1]*1.0);
676
677 // if the two check do not agree
678 if(check[0] != check[1]){
679 if(kfMother[0]) delete kfMother[0];
680 if(kfMother[1]) delete kfMother[1];
681 return kFALSE;
682 }
683
684 // now that the check[0] == check[1]
685 const Int_t type = check[0];
686
687 Float_t iMass =0.;
688 if(CheckSigns(v0)){
689 iMass = (type == 0) ? v0->GetEffMass(4, 2) : v0->GetEffMass(2, 4);
690 }
691 else{
692 iMass = (type == 0) ? v0->GetEffMass(2, 4) : v0->GetEffMass(4, 2);
693 }
694 Float_t iP = v0->P();
695
696 // Cuts
697 const Double_t cutCosPoint[2] = {0., 0.03}; // ORG [0., 0.03]
698 const Double_t cutDCA[2] = {0., 0.2}; // ORG [0., 0.2]
699 const Double_t cutProdVtxR[2] = {2., 30.}; // ORG [0., 24.]
700 const Double_t cutMass[2] = {1.11, 1.12}; // ORG [1.11, 1.12]
701 const Double_t cutChi2NDF = 5.; // ORG [5.]
702 // cundidate cuts
703 // opening angle as a function of L momentum
704 const Double_t cutOAngleP = 0.3 - 0.2*iP; // momentum dependent min. OAngle linear cut
705 // relative daughter momentum versusu mother momentum
706 // armenteros plot cuts
707 const Double_t cutQT = 0.03;
708 const Double_t cutAlpha = 0.7; // VERY strong - should supress the overlap with K0
709 // next cut see below
710
711 // compute the cut values
712
713 // cos pointing angle
714 Double_t cosPoint = v0->GetV0CosineOfPointingAngle();
715 cosPoint = TMath::ACos(cosPoint);
716
717 // DCA between daughters
718 Double_t dca = v0->GetDcaV0Daughters();
719
720 // Production vertex
721 Double_t x, y, z;
722 v0->GetXYZ(x,y,z);
723 Double_t r = TMath::Sqrt(x*x + y*y);
724
725 Int_t ix[2] = {0, 1};
726 if(1 == type){
727 ix[0] = 1;
728 ix[1] = 0;
729 }
730
731 // V0 chi2/ndf
732 Double_t chi2ndf = kfMother[type]->GetChi2()/kfMother[type]->GetNDF();
733
734 if(kfMother[0]) delete kfMother[0];
735 if(kfMother[1]) delete kfMother[1];
736
737 // Opening angle
738 Double_t oAngle = OpenAngle(v0);
739
740 // Relative daughter momentum
741 Double_t rP = (0 == check[0]) ? p[1]/p[0] : p[0]/p[1];
742
743 // Armenteros
744 Float_t ap[2];
745 Armenteros(v0, ap);
746
747 Double_t cutAP[2]; // a bit of woodoo :-)
748 cutAP[0] = 1.0 - (ap[0]-0.7 * ap[0]-0.7)*1.1 - 0.87;
749 cutAP[1] = 1.0 - (ap[0]+0.7 * ap[0]+0.7)*1.1 - 0.87;
750
751
752 //
753 // Apply the cuts, produce QA plots (with mass cut)
754 //
755
756 (type == 0) ? fQA->Fill("h_L_Mass", 0, iMass) : fQA->Fill("h_AL_Mass", 0, iMass);
757
758 if(iMass > cutMass[0] && iMass < cutMass[1]){
759 fQA->Fill("h_all_AP", 0, ap[0], ap[1]);
760 fQA->Fill("h_ProtonL_P", 0, p[ix[0]]);
761 fQA->Fill("h_PionL_P", 0, p[ix[1]]);
762 fQA->Fill("h_cut_L_Chi2", 0, chi2ndf);
763 fQA->Fill("h_cut_L_CosPoint", 0, cosPoint);
764 fQA->Fill("h_cut_L_CosPoint", 1, cosPoint);
765 fQA->Fill("h_cut_L_DCA", 0, dca);
766 fQA->Fill("h_cut_L_VtxR", 0, r);
767 fQA->Fill("h_cut_L_OAvP", 0, iP, oAngle);
768 fQA->Fill("h_cut_L_rdp_v_mp", 0, iP, rP);
769 if(0 ==type) fQA->Fill("h_cut_L_AP", 0, ap[0], ap[1]);
770 else fQA->Fill("h_cut_AL_AP", 0, ap[0], ap[1]);
771 fQA->Fill("h_cut_L_qT_v_mp", 0, iP, ap[1]);
772 }
773
774 if(cosPoint < cutCosPoint[0] || cosPoint > cutCosPoint[1]) return kFALSE;
775 (type == 0) ? fQA->Fill("h_L_Mass", 1, iMass) : fQA->Fill("h_AL_Mass", 1, iMass);
776 if(iMass > cutMass[0] && iMass < cutMass[1]){
777 fQA->Fill("h_ProtonL_P", 1, p[ix[0]]);
778 fQA->Fill("h_PionL_P", 1, p[ix[1]]);
779 fQA->Fill("h_cut_L_DCA", 1, dca);
780 }
781
782 if(dca < cutDCA[0] || dca > cutDCA[1]) return kFALSE;
783 (type == 0) ? fQA->Fill("h_L_Mass", 2, iMass) : fQA->Fill("h_AL_Mass", 2, iMass);
784 if(iMass > cutMass[0] && iMass < cutMass[1]){
785 fQA->Fill("h_ProtonL_P", 2, p[ix[0]]);
786 fQA->Fill("h_PionL_P", 2, p[ix[1]]);
787 fQA->Fill("h_cut_L_VtxR", 1, r);
788 }
789
790 if(r < cutProdVtxR[0] || r > cutProdVtxR[1]) return kFALSE;
791 (type == 0) ? fQA->Fill("h_L_Mass", 3, iMass) : fQA->Fill("h_AL_Mass", 3, iMass);
792 if(iMass > cutMass[0] && iMass < cutMass[1]){
793 fQA->Fill("h_ProtonL_P", 3, p[ix[0]]);
794 fQA->Fill("h_PionL_P", 3, p[ix[1]]);
795 fQA->Fill("h_cut_L_Chi2", 1, chi2ndf);
796 }
797
798
799 if(chi2ndf > cutChi2NDF) return kFALSE;
800 (type == 0) ? fQA->Fill("h_L_Mass", 4, iMass) : fQA->Fill("h_AL_Mass", 4, iMass);
801 if(iMass > cutMass[0] && iMass < cutMass[1]){
802 fQA->Fill("h_ProtonL_P", 4, p[ix[0]]);
803 fQA->Fill("h_PionL_P", 4, p[ix[1]]);
804 fQA->Fill("h_cut_L_OAvP", 1, iP, oAngle);
805 }
806
807 if(oAngle < cutOAngleP) return kFALSE;
808 (type == 0) ? fQA->Fill("h_L_Mass", 5, iMass) : fQA->Fill("h_AL_Mass", 5, iMass);
809 if(iMass > cutMass[0] && iMass < cutMass[1]){
810 fQA->Fill("h_ProtonL_P", 5, p[ix[0]]);
811 fQA->Fill("h_PionL_P", 5, p[ix[1]]);
812 fQA->Fill("h_cut_L_rdp_v_mp", 1, iP, rP);
813 if(0 == type) fQA->Fill("h_cut_L_AP", 1, ap[0], ap[1]);
814 else fQA->Fill("h_cut_AL_AP", 1, ap[0], ap[1]);
815 fQA->Fill("h_cut_L_qT_v_mp", 1, iP, ap[1]);
816 fQA->Fill("h_all_AP", 1, ap[0], ap[1]);
817 }
818
819 if(ap[1] < cutQT) return kFALSE;
820 (type == 0) ? fQA->Fill("h_L_Mass", 6, iMass) : fQA->Fill("h_AL_Mass", 6, iMass);
821 if(iMass > cutMass[0] && iMass < cutMass[1]){
822 fQA->Fill("h_ProtonL_P", 6, p[ix[0]]);
823 fQA->Fill("h_PionL_P", 6, p[ix[1]]);
824 }
825
826 if(ap[1] > cutAP[type]) return kFALSE;
827 (type == 0) ? fQA->Fill("h_L_Mass", 7, iMass) : fQA->Fill("h_AL_Mass", 7, iMass);
828 if(iMass > cutMass[0] && iMass < cutMass[1]){
829 fQA->Fill("h_ProtonL_P", 7, p[ix[0]]);
830 fQA->Fill("h_PionL_P", 7, p[ix[1]]);
831 }
832 if(TMath::Abs(ap[0]) > cutAlpha) return kFALSE;
833 (type == 0) ? fQA->Fill("h_L_Mass", 8, iMass) : fQA->Fill("h_AL_Mass", 8, iMass);
834 if(iMass > cutMass[0] && iMass < cutMass[1]){
835 fQA->Fill("h_ProtonL_P", 8, p[ix[0]]);
836 fQA->Fill("h_PionL_P", 8, p[ix[1]]);
837 }
838
839 if(iMass < cutMass[0] || iMass > cutMass[1]) {
840 return kFALSE;
841 }
842
843 // all cuts passed
844
845 // assign the lambda type value: Lambda: kTRUE, Anti-Lambda: kFALSE
846 isLambda = (0 == type) ? kTRUE : kFALSE;
847
848
849 return kTRUE;
850}
851//________________________________________________________________
852Double_t AliHFEV0cuts::OpenAngle(AliESDv0 *v0) const {
853 //
854 // Opening angle between two daughter tracks
855 //
856 Double_t mn[3] = {0,0,0};
857 Double_t mp[3] = {0,0,0};
858
859
860 v0->GetNPxPyPz(mn[0],mn[1],mn[2]);//reconstructed cartesian momentum components of negative daughter;
861 v0->GetPPxPyPz(mp[0],mp[1],mp[2]);//reconstructed cartesian momentum components of positive daughter;
862
863
864 Double_t openAngle = 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])));
865
866 return TMath::Abs(openAngle);
867}
868//________________________________________________________________
869Double_t AliHFEV0cuts::PsiPair(AliESDv0 *v0) {
870 //
871 // Angle between daughter momentum plane and plane
872 //
873
874 if(!fInputEvent) return -1.;
875
876 Float_t magField = fInputEvent->GetMagneticField();
877
878 Int_t pIndex = -1;
879 Int_t nIndex = -1;
880 if(CheckSigns(v0)){
881 pIndex = v0->GetPindex();
882 nIndex = v0->GetNindex();
883 }
884 else{
885 pIndex = v0->GetNindex();
886 nIndex = v0->GetPindex();
887 }
888
889
890 AliESDtrack* daughter[2];
891
892 daughter[0] = dynamic_cast<AliESDtrack *>(fInputEvent->GetTrack(pIndex));
893 daughter[1] = dynamic_cast<AliESDtrack *>(fInputEvent->GetTrack(nIndex));
894
895 Double_t x, y, z;
896 v0->GetXYZ(x,y,z);//Reconstructed coordinates of V0; to be replaced by Markus Rammler's method in case of conversions!
897
898 Double_t mn[3] = {0,0,0};
899 Double_t mp[3] = {0,0,0};
900
901
902 v0->GetNPxPyPz(mn[0],mn[1],mn[2]);//reconstructed cartesian momentum components of negative daughter;
903 v0->GetPPxPyPz(mp[0],mp[1],mp[2]);//reconstructed cartesian momentum components of positive daughter;
904
905
906 Double_t deltat = 1.;
907 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
908
909 Double_t radiussum = TMath::Sqrt(x*x + y*y) + 50;//radius to which tracks shall be propagated
910
911 Double_t momPosProp[3];
912 Double_t momNegProp[3];
913
914 AliExternalTrackParam pt(*daughter[0]), nt(*daughter[1]);
915
916 Double_t psiPair = 4.;
917
918 if(nt.PropagateTo(radiussum,magField) == 0)//propagate tracks to the outside
919 psiPair = -5.;
920 if(pt.PropagateTo(radiussum,magField) == 0)
921 psiPair = -5.;
922 pt.GetPxPyPz(momPosProp);//Get momentum vectors of tracks after propagation
923 nt.GetPxPyPz(momNegProp);
924
925 Double_t pEle =
926 TMath::Sqrt(momNegProp[0]*momNegProp[0]+momNegProp[1]*momNegProp[1]+momNegProp[2]*momNegProp[2]);//absolute momentum value of negative daughter
927 Double_t pPos =
928 TMath::Sqrt(momPosProp[0]*momPosProp[0]+momPosProp[1]*momPosProp[1]+momPosProp[2]*momPosProp[2]);//absolute momentum value of positive daughter
929
930 Double_t scalarproduct =
931 momPosProp[0]*momNegProp[0]+momPosProp[1]*momNegProp[1]+momPosProp[2]*momNegProp[2];//scalar product of propagated positive and negative daughters' momenta
932
933 Double_t chipair = TMath::ACos(scalarproduct/(pEle*pPos));//Angle between propagated daughter tracks
934
935 psiPair = TMath::Abs(TMath::ASin(deltat/chipair));
936
937 return psiPair;
938}
939//________________________________________________________________
940AliKFParticle *AliHFEV0cuts::CreateMotherParticle(AliVTrack* const pdaughter, AliVTrack* const ndaughter, Int_t pspec, Int_t nspec){
941 //
942 // Creates a mother particle
943 //
944 AliKFParticle pkfdaughter(*pdaughter, pspec);
945 AliKFParticle nkfdaughter(*ndaughter, nspec);
946
947 // - check if the daughter particles are coming from the primary vertex
948 // - check the number of tracks that take part in the creaton of primary vertex.
949 // important: after removeal of candidate tracks there must be at least 2 tracks left
950 // otherwise the primary vertex will be corrupted
951
952 // ESD Analyis
953 //const AliESDVertex *esdvertex = dynamic_cast<const AliESDVertex *>(fInputEvent->GetPrimaryVertex());
954 //if(!esdvertex) return NULL;
955 //UShort_t *contrib = esdvertex->GetIndices();
956
957 //
958 // not using the removal of the daughter track now
959 //
960 // Int_t nTracks = esdvertex->GetNIndices();
961 // printf(" -D: N Vertex tracks: %i\n", nTracks);
962 // printf(" -D: N Contributors: %i\n", fPrimaryVertex->GetNContributors());
963 // Int_t nfound = 0;
964 // for(Int_t id = 0; id < esdvertex->GetNIndices(); id++){
965 // if(contrib[id] == pdaughter->GetID()){
966 // if( (nTracks - nfound) <= 2 ) return NULL;
967 // *fPrimaryVertex -= pkfdaughter;
968 // removed[0] = kTRUE;
969 // nfound++;
970 // }
971 // if(contrib[id] == ndaughter->GetID()){
972 // if( (nTracks - nfound) <=2 ) return NULL;
973 // *fPrimaryVertex -= nkfdaughter;
974 // removed[1] = kTRUE;
975 // nfound++;
976 // }
977 // if(nfound == 2) break;
978 // }
979
980 // printf(" -D: n removed: %i\n", nfound);
981
982 // Create the mother particle
983 AliKFParticle *m = new AliKFParticle(pkfdaughter, nkfdaughter);
984
985 AliKFVertex improvedVertex = *fPrimaryVertex;
986 improvedVertex += *m;
987 m->SetProductionVertex(improvedVertex);
988
989 // update 15/06/2010
990 // mother particle will not be added to primary vertex but only to its copy
991 // as this confilcts with calling
992 // m->SetPrimaryVertex() function and
993 // subsequently removing the mother particle afterwards
994 // Sourse: Sergey Gorbunov
995
996 return m;
997}
998//_________________________________________________
999Bool_t AliHFEV0cuts::LooseRejectK0(AliESDv0 * const v0) const {
1000 //
1001 // Reject K0 based on loose cuts
1002 //
1003 Double_t mass = v0->GetEffMass(AliPID::kPion, AliPID::kPion);
1004 if(mass > 0.494 && mass < 0.501) return kTRUE;
1005 return kFALSE;
1006}
1007
1008//_________________________________________________
1009Bool_t AliHFEV0cuts::LooseRejectLambda(AliESDv0 * const v0) const {
1010 //
1011 // Reject Lambda based on loose cuts
1012 //
1013 Double_t mass1 = v0->GetEffMass(AliPID::kPion, AliPID::kProton);
1014 Double_t mass2 = v0->GetEffMass(AliPID::kProton, AliPID::kPion);
1015
1016 if(mass1 > 1.1 && mass1 < 1.12) return kTRUE;
1017 if(mass2 > 1.1 && mass2 < 1.12) return kTRUE;
1018 return kFALSE;
1019}
1020
1021//_________________________________________________
1022Bool_t AliHFEV0cuts::LooseRejectGamma(AliESDv0 * const v0) const {
1023 //
1024 // Reject Lambda based on loose cuts
1025 //
1026
1027 Double_t mass = v0->GetEffMass(AliPID::kElectron, AliPID::kElectron);
1028
1029 if(mass < 0.02) return kTRUE;
1030 return kFALSE;
1031}
1032//___________________________________________________________________
1033void AliHFEV0cuts::Armenteros(AliESDv0 *v0, Float_t val[2]){
1034 //
1035 // computes the Armenteros variables for given V0
1036 // fills the histogram
1037 // returns the values via "val"
1038 //
1039
1040 Double_t mn[3] = {0,0,0};
1041 Double_t mp[3] = {0,0,0};
1042 Double_t mm[3] = {0,0,0};
1043
1044 if(CheckSigns(v0)){
1045 v0->GetNPxPyPz(mn[0],mn[1],mn[2]); //reconstructed cartesian momentum components of negative daughter
1046 v0->GetPPxPyPz(mp[0],mp[1],mp[2]); //reconstructed cartesian momentum components of positive daughter
1047 }
1048 else{
1049 v0->GetPPxPyPz(mn[0],mn[1],mn[2]); //reconstructed cartesian momentum components of negative daughter
1050 v0->GetNPxPyPz(mp[0],mp[1],mp[2]); //reconstructed cartesian momentum components of positive daughter
1051 }
1052 v0->GetPxPyPz(mm[0],mm[1],mm[2]); //reconstructed cartesian momentum components of mother
1053
1054 TVector3 vecN(mn[0],mn[1],mn[2]);
1055 TVector3 vecP(mp[0],mp[1],mp[2]);
1056 TVector3 vecM(mm[0],mm[1],mm[2]);
1057
1058 Double_t thetaP = acos((vecP * vecM)/(vecP.Mag() * vecM.Mag()));
1059 Double_t thetaN = acos((vecN * vecM)/(vecN.Mag() * vecM.Mag()));
1060
1061 Double_t alfa = ((vecP.Mag())*cos(thetaP)-(vecN.Mag())*cos(thetaN))/
1062 ((vecP.Mag())*cos(thetaP)+(vecN.Mag())*cos(thetaN)) ;
1063 Double_t qt = vecP.Mag()*sin(thetaP);
1064
1065 val[0] = alfa;
1066 val[1] = qt;
1067
1068}
1069//___________________________________________________________________
1070Bool_t AliHFEV0cuts::CheckSigns(AliESDv0* const v0){
1071 //
1072 // check wheter the sign was correctly applied to
1073 // V0 daughter tracks
1074 //
1075
1076 Bool_t correct = kFALSE;
1077
1078 Int_t pIndex = 0, nIndex = 0;
1079 pIndex = v0->GetPindex();
1080 nIndex = v0->GetNindex();
1081
1082 AliESDtrack* d[2];
1083 d[0] = dynamic_cast<AliESDtrack*>(fInputEvent->GetTrack(pIndex));
1084 d[1] = dynamic_cast<AliESDtrack*>(fInputEvent->GetTrack(nIndex));
1085
1086 Int_t sign[2];
1087 sign[0] = (int)d[0]->GetSign();
1088 sign[1] = (int)d[1]->GetSign();
1089
1090 if(-1 == sign[0] && 1 == sign[1]){
1091 correct = kFALSE;
1092 //v0->SetIndex(0, pIndex); // set the index of the negative v0 track
1093 //v0->SetIndex(1, nIndex); // set the index of the positive v0 track
1094 }
1095 else{
1096 correct = kTRUE;
1097 }
1098
1099 //pIndex = v0->GetPindex();
1100 //nIndex = v0->GetNindex();
1101 //printf("-D2: P: %i, N: %i\n", pIndex, nIndex);
1102
1103 return correct;
1104}