]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/KINK/AliAnalysisKinkESDMC.cxx
Particle mass access updated: Evi Ganoti, University of Athens (pganoti@phys.uoa.gr)
[u/mrichter/AliRoot.git] / PWG2 / KINK / AliAnalysisKinkESDMC.cxx
1 /**************************************************************************
2  * Authors: Martha Spyropoulou-Stassinaki and the  members 
3  * of the Greek group at Physics Department of Athens University
4  * Paraskevi Ganoti, Anastasia Belogianni and Filimon Roukoutakis 
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 //-----------------------------------------------------------------
17 //                 AliAnalysisKinkESDMC class
18 //       Example of an analysis task for kink topology study
19 //      Kaons from kink topology are 'identified' in this code
20 //-----------------------------------------------------------------
21
22 #include "TChain.h"
23 #include "TTree.h"
24 #include "TH1F.h"
25 #include "TH2F.h"
26 #include "TH1D.h"
27 #include "TH2D.h"
28 #include "TParticle.h"
29 #include <TVector3.h>
30 #include "TF1.h"
31
32 #include "AliAnalysisTaskSE.h"
33 #include "AliAnalysisManager.h"
34
35 #include "AliESDEvent.h"
36 #include "AliESDInputHandler.h"
37
38 #include "AliMCEventHandler.h"
39 #include "AliMCEvent.h"
40
41 #include "AliAnalysisKinkESDMC.h"
42
43 #include "AliStack.h"
44 #include "AliESDpid.h"
45 #include "AliESDkink.h"
46
47 ClassImp(AliAnalysisKinkESDMC)
48
49 //________________________________________________________________________
50 AliAnalysisKinkESDMC::AliAnalysisKinkESDMC() 
51   : AliAnalysisTaskSE(), fESD(0),fHistPtESD(),fHistPt(0),fHistQtAll(0),fHistQt1(0),fHistQt2(0)
52   , fHistPtKaon(0),fHistPtKPDG(0),fHistEta(0),fHistEtaK(0),fptKMC(0),fMultiplMC(0),fESDMult(0),fgenpt(0),frad(0),
53   fKinkKaon(0), fKinkKaonBg(0), fM1kaon(0),  fgenPtEtR(0),fPtKink(0),  fptKink(0),
54    fcodeH(0), fdcodeH(0), fAngMomK(0),fAngMomPi(0),
55   fRpr(0),fZpr(0),
56    fAngMomKC(0), 
57   fZvXv(0),fZvYv(0), fXvYv(0),  fPtPrKink(0), f1(0),  f2(0),
58       fListOfHistos(0)
59
60 {
61   // Constructor
62
63 }
64
65
66 //________________________________________________________________________
67 AliAnalysisKinkESDMC::AliAnalysisKinkESDMC(const char *name) 
68   : AliAnalysisTaskSE(name), fESD(0),fHistPtESD(0),fHistPt(0),fHistQtAll(0),fHistQt1(0),fHistQt2(0)
69   , fHistPtKaon(0),fHistPtKPDG(0),fHistEta(0),fHistEtaK(0),fptKMC(0),fMultiplMC(0),fESDMult(0),fgenpt(0),frad(0),
70   fKinkKaon(0), fKinkKaonBg(0), fM1kaon(0),  fgenPtEtR(0),fPtKink(0),  fptKink(0),
71    fcodeH(0), fdcodeH(0), fAngMomK(0),fAngMomPi(0),
72   fRpr(0),fZpr(0),
73    fAngMomKC(0),
74    fZvXv(0),fZvYv(0), fXvYv(0), fPtPrKink(0),  f1(0), f2(0),
75       fListOfHistos(0)
76   // , faod(0), fSignalTree(0)
77
78 {
79   // Constructor
80
81   // Define input and output slots here
82   // Input slot #0 works with a TChain
83   DefineInput(0, TChain::Class());
84   DefineOutput(1, TList::Class());
85 }
86
87 //________________________________________________________________________
88 void AliAnalysisKinkESDMC::ConnectInputData(Option_t *) 
89 {
90   // Connect ESD or AOD here
91   // Called once
92
93   TTree* tree = dynamic_cast<TTree*> (GetInputData(0));
94   if (!tree) {
95     Printf("ERROR: Could not read chain from input slot 0");
96   } else {
97
98     AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
99
100     if (!esdH) {
101       Printf("ERROR: Could not get ESDInputHandler");
102     } else
103       fESD = esdH->GetEvent();
104   }
105 }
106
107 //________________________________________________________________________
108 void AliAnalysisKinkESDMC::CreateOutputObjects() 
109 {
110   // Create histograms
111   // Called once
112   
113    f1=new TF1("f1","((atan([0]*[1]*(1.0/(sqrt((x^2)*(1.0-([1]^2))-([0]^2)*([1]^2))))))*180.)/[2]",1.1,10.0);
114    f1->SetParameter(0,0.493677);
115    f1->SetParameter(1,0.9127037);
116    f1->SetParameter(2,TMath::Pi());
117
118
119    f2=new TF1("f2","((atan([0]*[1]*(1.0/(sqrt((x^2)*(1.0-([1]^2))-([0]^2)*([1]^2))))))*180.)/[2]",0.1,10.0);
120    f2->SetParameter(0,0.13957018);
121    f2->SetParameter(1,0.2731374);
122    f2->SetParameter(2,TMath::Pi());
123    //Open file  1= CAF 
124     //OpenFile(1); 
125
126   fHistPtESD = new TH1F("fHistPtESD", "P_{T} distribution", 60, 0.0, 6.0);
127   fHistPtESD->GetXaxis()->SetTitle("P_{T} (GeV/c)");
128   fHistPtESD->GetYaxis()->SetTitle("dN/dP_{T} (c/GeV)");
129   fHistPtESD->SetMarkerStyle(kFullCircle);
130   fHistPt = new TH1F("fHistPt", "P_{T} distribution", 60, 0.0, 6.0); 
131   fHistQtAll = new TH1F("fHistQtAll", "Q_{T} distr All Kinks ",100, 0.0,.250); 
132   fHistQt1= new TH1F("fHistQt1", "Q_{T} distribution",100, 0.0,.300); 
133   fHistQt2= new TH1F("fHistQt2", "Q_{T} distribution",100, 0.0,.300); 
134   fHistPtKaon = new TH1F("fHistPtKaon", "P_{T}Kaon distribution", 60, 0.0, 6.0); 
135   fHistPtKPDG = new TH1F("fHistPtKPDG", "P_{T}Kaon distribution", 60, 0.0, 6.0); 
136   fHistEta= new TH1F("fHistEta", "Eta distribution", 26,-1.3, 1.3); 
137   fHistEtaK= new TH1F("fHistEtaK", "EtaK distribution", 26,-1.3, 1.3); 
138   fptKMC= new TH1F("fptKMC", "P_{T}Kaon generated", 60, 0.0, 6.0); 
139   fMultiplMC= new TH1F("fMultiplMC", "charge multiplicity MC", 60, 0.5,120.5); 
140   fESDMult= new TH1F("fESDMult", "charge multipliESD", 60, 0.5,120.5); 
141   fgenpt= new TH1F("fgenpt", "genpt   K distribution", 60, 0.0, 6.0); 
142   frad= new TH1F("frad", "radius  K generated",100,0.,400.); 
143   fKinkKaon= new TH1F("fKinkKaon", "P_{T}Kaon kinks identi", 60, 0.0, 6.0); 
144   fKinkKaonBg= new TH1F("fKinkKaonBg", "P_{T}Kaon kinks backgr", 60, 0.0, 6.0); 
145   fM1kaon= new TH1F("fM1kaon","Invar m(kaon) from kink->mu+netrino decay",80,0.0, 0.8); 
146   fgenPtEtR= new TH1F("fgenPtEtR", "P_{T}Kaon distribution", 60, 0.0, 6.0); 
147   fPtKink= new TH1F("fPtKink", "P_{T}Kaon Kink  bution", 60, 0.0, 6.0); 
148   fptKink= new TH1F("fptKink", "P_{T}Kaon Kink  bution", 60, 0.0, 6.0); 
149   fcodeH   = new TH2F("fcodeH", "code vrs dcode dist. kinks,K",100,0.,2500.,100,0.,2500.);
150   fdcodeH = new TH2F("fdcodeH", "code vrs dcode dist. kinks,K",100,0.,2500.,100,0.,2500.);
151   fAngMomK= new TH2F("fAngMomK","Decay angle vrs Mother Mom,K",100,0.0,5.0,80,0.,80.);
152   fAngMomPi= new TH2F("fAngMomPi","Decay angle vrs Mother Mom,Pi",100,0.0,5.0,80,0.,80.);
153
154   fRpr = new TH1D("fRpr", "rad distribution  PID pr",50,0.0, 2.5);
155   fZpr = new TH1D("fZpr", "z distribution PID pr  ",60,-15.,15.);
156   fAngMomKC= new TH2F("fAngMomKC","Decay angle vrs Mother Mom,K",100,0.0,5.0,80,0.,80.);
157   fZvXv= new TH2F("fZvXv","Xv-Zv main vtx",60,-0.5,0.5,60, -15., 15.0);
158   fZvYv= new TH2F("fZvYv","Yv-Zv main vtx",60,-0.5,0.5, 60, -15., 15.);
159   fXvYv= new TH2F("fXvYv","Xv-Yv main vtx", 60,-1.5,1.5, 60, -1.5, 1.5);
160   fPtPrKink=new TH1F("fPtPrKink","pt of ESD  kaonKink tracks", 60, 0.0, 6.0);
161
162    fListOfHistos=new TList();
163
164    fListOfHistos->Add(fHistPtESD);
165    fListOfHistos->Add(fHistPt);
166    fListOfHistos->Add(fHistQtAll);
167    fListOfHistos->Add(fHistQt1);
168    fListOfHistos->Add(fHistQt2);
169    fListOfHistos->Add(fHistPtKaon);
170    fListOfHistos->Add(fHistPtKPDG);
171    fListOfHistos->Add(fHistEta);
172    fListOfHistos->Add(fHistEtaK);
173    fListOfHistos->Add(fptKMC);
174    fListOfHistos->Add(fMultiplMC);
175    fListOfHistos->Add(fESDMult);
176    fListOfHistos->Add(fgenpt);
177    fListOfHistos->Add(frad);
178    fListOfHistos->Add(fKinkKaon);
179    fListOfHistos->Add(fKinkKaonBg);
180    fListOfHistos->Add(fM1kaon);
181    fListOfHistos->Add(fgenPtEtR);
182    fListOfHistos->Add(fPtKink);
183    fListOfHistos->Add(fptKink);
184    fListOfHistos->Add(fcodeH);
185    fListOfHistos->Add(fdcodeH);
186    fListOfHistos->Add(fdcodeH);
187    fListOfHistos->Add(fAngMomK);
188    fListOfHistos->Add(fAngMomPi);
189    fListOfHistos->Add(fRpr);
190    fListOfHistos->Add(fZpr);
191    fListOfHistos->Add(fAngMomKC);
192    fListOfHistos->Add(fZvXv);
193    fListOfHistos->Add(fZvYv);
194    fListOfHistos->Add(fXvYv);
195    fListOfHistos->Add(fPtPrKink);
196 }
197
198 //________________________________________________________________________
199 void AliAnalysisKinkESDMC::Exec(Option_t *) 
200 {
201   // Main loop
202   // Called for each event
203
204   // Process MC truth, therefore we receive the AliAnalysisManager and ask it for the AliMCEventHandler
205   // This handler can return the current MC event
206  
207   AliMCEventHandler* eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
208   if (!eventHandler) {
209     Printf("ERROR: Could not retrieve MC event handler");
210     return;
211   }
212
213   AliMCEvent* mcEvent = eventHandler->MCEvent();
214   if (!mcEvent) {
215      Printf("ERROR: Could not retrieve MC event");
216      return;
217   }
218
219   Printf("MC particles: %d", mcEvent->GetNumberOfTracks());
220
221   if (!fESD) {
222     Printf("ERROR: fESD not available");
223     return;
224   }
225   AliStack* stack=mcEvent->Stack();
226   
227 //primary tracks  in MC
228          Int_t  nPrim = stack->GetNprimary();
229 //
230   const AliESDVertex *vertex=GetEventVertex(fESD);
231     if(!vertex) return;
232 //
233   Double_t vpos[3];
234   vertex->GetXYZ(vpos);
235     fZpr->Fill(vpos[2]);         
236  fZvXv->Fill( vpos[1], vpos[2]);  
237  fZvYv->Fill( vpos[0], vpos[2]);  
238
239   Double_t vtrack[3], ptrack[3];
240   
241      
242
243    Int_t nGoodTracks =  fESD->GetNumberOfTracks();
244     fESDMult->Fill(nGoodTracks);
245 //
246 // track loop
247   for (Int_t iTracks = 0; iTracks < fESD->GetNumberOfTracks(); iTracks++) {
248     AliESDtrack* track = fESD->GetTrack(iTracks);
249     if (!track) {
250       Printf("ERROR: Could not receive track %d", iTracks);
251       continue;
252     }
253     
254
255     fHistPt->Fill(track->Pt());
256
257       //   pt cut at 300 MeV
258         if ( (track->Pt())<.300)continue;
259
260     UInt_t status=track->GetStatus();
261
262   //  if((status&AliESDtrack::kITSrefit)==0) continue;
263     if((status&AliESDtrack::kTPCrefit)==0) continue;
264       if((track->GetTPCchi2()/track->GetTPCclusters(0))>3.5) continue;
265
266       Double_t extCovPos[15];
267       track->GetExternalCovariance(extCovPos);    
268       if(extCovPos[0]>2) continue;
269       if(extCovPos[2]>2) continue;    
270       if(extCovPos[5]>0.5) continue;  
271       if(extCovPos[9]>0.5) continue;
272       if(extCovPos[14]>2) continue;
273
274
275     track->GetXYZ(vtrack);
276  fXvYv->Fill(vtrack[0],vtrack[1]);  
277
278 // track momentum
279      track->GetPxPyPz(ptrack);
280     
281     TVector3 trackMom(ptrack[0],ptrack[1],ptrack[2]);
282     
283     Double_t trackEta=trackMom.Eta();
284     
285     
286     
287     Float_t nSigmaToVertex = GetSigmaToVertex(track);  
288       
289     Float_t bpos[2];
290     Float_t bCovpos[3];
291     track->GetImpactParameters(bpos,bCovpos);
292     
293     if (bCovpos[0]<=0 || bCovpos[2]<=0) {
294      Printf("Estimated b resolution lower or equal zero!");
295      bCovpos[0]=0; bCovpos[2]=0;
296     }
297
298     Float_t dcaToVertexXYpos = bpos[0];
299     Float_t dcaToVertexZpos = bpos[1];
300     
301     fRpr->Fill(dcaToVertexZpos);
302
303     if(nSigmaToVertex>=4) continue;
304     if((dcaToVertexXYpos>3.0)||(dcaToVertexZpos>3.0)) continue;
305     
306     
307     Int_t label = track->GetLabel();
308     label = TMath::Abs(label);
309     TParticle * part = stack->Particle(label);
310     if (!part) continue;
311
312  //  cut on eta 
313         if(  (TMath::Abs(trackEta )) > 1.1 ) continue;
314     fHistPtESD->Fill(track->Pt());
315
316    // Add Kink analysis
317    
318                 Int_t indexKinkPos=track->GetKinkIndex(0);
319 //  loop on kinks
320                 if(indexKinkPos<0){
321                fPtKink->Fill(track->Pt()); /// pt from track
322
323         // select kink class    
324
325           AliESDkink *kink=fESD->GetKink(TMath::Abs(indexKinkPos)-1);
326 //
327         
328           Int_t eSDfLabel1=kink->GetLabel(0);
329           TParticle *particle1= stack->Particle(TMath::Abs(eSDfLabel1));
330           Int_t code1= particle1->GetPdgCode();
331           
332           Int_t eSDfLabeld=kink->GetLabel(1);
333           TParticle *particled= stack->Particle(TMath::Abs(eSDfLabeld));
334           Int_t dcode1= particled->GetPdgCode();
335
336            const TVector3 motherMfromKink(kink->GetMotherP());
337            const TVector3 daughterMKink(kink->GetDaughterP());
338            Float_t qT=kink->GetQt();
339
340            fHistEta->Fill(trackEta) ;  //   Eta distr
341            fHistQtAll->Fill(qT) ;  //  Qt   distr
342                   
343            fptKink->Fill(motherMfromKink.Pt()); /// pt from kink
344
345 //       fake kinks, small Qt and small kink angle
346            if ( qT<0.012) fcodeH->Fill( TMath::Abs(code1), TMath::Abs(dcode1) );
347            if(qT<0.012) continue;  // remove fake kinks
348
349            Float_t kinkAngle=TMath::RadToDeg()*kink->GetAngle(2);
350          
351 //remove the double taracks 
352             
353            if( (kinkAngle<1.)  ) continue;
354 //
355            if( (qT>0.05)&& ( qT<0.25)  )   fHistQt2->Fill(qT);  // candidate kaon kinks
356
357            if (TMath::Abs(code1)==321)     fHistPtKPDG->Fill(track->Pt());  // ALL KAONS (pdg) inside kink sample
358
359 //          maximum decay angle at a given mother momentum
360            Double_t maxDecAngKmu=f1->Eval(motherMfromKink.Mag(),0.,0.,0.);
361            Double_t maxDecAngpimu=f2->Eval(motherMfromKink.Mag(),0.,0.,0.);
362 //  two dimensional plot 
363                 if(TMath::Abs(code1)==321) fAngMomK->Fill(motherMfromKink.Mag(), kinkAngle); 
364                 if(TMath::Abs(code1)==211) fAngMomPi->Fill(motherMfromKink.Mag(), kinkAngle); 
365 //
366 // invariant mass of mother track decaying to mu
367          Float_t energyDaughterMu=TMath::Sqrt(daughterMKink.Mag()*daughterMKink.Mag()+0.105658*0.105658);
368          Float_t p1XM= motherMfromKink.Px();
369          Float_t p1YM= motherMfromKink.Py();
370          Float_t p1ZM= motherMfromKink.Pz();
371          Float_t p2XM= daughterMKink.Px();
372          Float_t p2YM= daughterMKink.Py();
373          Float_t p2ZM= daughterMKink.Pz();
374          Float_t p3Daughter=TMath::Sqrt(((p1XM-p2XM)*(p1XM-p2XM))+((p1YM-p2YM)*(p1YM-p2YM))+((p1ZM-p2ZM)*(p1ZM-p2ZM)));
375          Double_t invariantMassKmu= TMath::Sqrt((energyDaughterMu+p3Daughter)*(energyDaughterMu+p3Daughter)-motherMfromKink.Mag()*motherMfromKink.Mag());
376   
377          fM1kaon->Fill(invariantMassKmu);
378
379 //  kaon selection from kinks
380 if((kinkAngle>maxDecAngpimu)&&(qT>0.05)&&(qT<0.25)&&((kink->GetR()>110.)&&(kink->GetR()<230.))&&(TMath::Abs(trackEta)<1.1)&&(invariantMassKmu<0.6)) {
381                 fHistEtaK->Fill(trackEta);
382
383                 if(TMath::Abs(code1)==321) fAngMomKC->Fill(motherMfromKink.Mag(), kinkAngle);  // real kaons
384
385                 fHistQt1  ->Fill(qT) ;  //  Qt   distr
386                 if ( (kinkAngle>maxDecAngKmu) && ( motherMfromKink.Mag()> 1.1 ) ) continue; // maximum angle selection
387
388                                     fHistPtKaon->Fill(track->Pt());   //all PID kink-kaon
389
390 // background inside the identified kaons, e.g KK  ,  Kp
391          if((TMath::Abs(code1)==321)&&(( TMath::Abs(dcode1)) >=(TMath::Abs(code1)) )) fdcodeH->Fill( TMath::Abs(code1), TMath::Abs(dcode1));
392
393          if((TMath::Abs(code1)==321) ) { 
394
395          if ( label< nPrim ) fPtPrKink->Fill(track->Pt());
396              fKinkKaon->Fill(track->Pt());        }
397          else {
398              fKinkKaonBg->Fill(track->Pt());     
399           }   // primary and secondary 
400
401         }  //  kink selection 
402                   
403
404         }  //End Kink Information    
405   
406
407   } //track loop 
408
409   // loop over mc particles
410
411   // variable to count tracks
412   Int_t nMCTracks = 0;
413
414 for (Int_t iMc = 0; iMc < nPrim; ++iMc)
415   {
416    // AliDebug(AliLog::kDebug+1, Form("MC Loop: Processing particle %d.", iMc));
417
418     TParticle* particle = stack->Particle(iMc);
419
420     if (!particle)
421     {
422     //  AliDebug(AliLog::kError, Form("UNEXPECTED: particle with label %d not found in stack (mc loop).", iMc));
423       continue;
424     }
425
426  //   
427     Float_t eta = particle->Eta();
428
429     if (TMath::Abs(eta) > 1.1)
430       continue;
431
432             Double_t ptK = particle->Pt();
433
434    if( ptK <0.300) continue;
435
436       Float_t code = particle->GetPdgCode();
437       if (code==321 || code==-321){
438             
439               fptKMC->Fill(ptK);
440     
441
442             Int_t firstD=particle->GetFirstDaughter();
443             Int_t lastD=particle->GetLastDaughter();
444 //loop on secondaries
445              for (Int_t k=firstD;k<lastD;k++) {
446              TParticle*    daughter1=stack->Particle(k+1);
447              Float_t dcode = daughter1->GetPdgCode();
448
449     
450
451                  frad->Fill(daughter1->R());
452
453                  if (((daughter1->R())>110)&&((daughter1->R())<230) ){
454         if (( ( code==321 )&& ( dcode ==-13  ))|| (( code == -321 )&& ( dcode == 13))) fgenPtEtR->Fill( ptK );//to muon
455         if (( ( code==321 )&& ( dcode ==211  ))|| (( code == -321 )&& ( dcode ==-211))) fgenPtEtR->Fill( ptK );//to pion
456      //   if (( ( code==321 )&& ( dcode ==-11  ))|| (( code == -321 )&& ( dcode ==11))) fgenPtEtR->Fill( ptK );// to electr
457
458        }
459 //  next only  to mu decay
460                  if (((code==321)&&(dcode==-13))||((code==-321)&&(dcode==13))){
461
462                  
463                     if (((daughter1->R())>110)&&((daughter1->R())<230)){
464                     fgenpt->Fill(ptK);
465                     }
466                     }
467        }
468
469       }
470
471
472       
473     nMCTracks++;
474   }// end of mc particle
475
476 //  printf("hello mc");
477   fMultiplMC->Fill(nMCTracks);
478
479   PostData(1, fListOfHistos);
480
481 }      
482
483 //________________________________________________________________________
484 void AliAnalysisKinkESDMC::Terminate(Option_t *) 
485 {
486   // Draw result to the screen 
487   // Called once at the end of the query
488
489 }
490 //____________________________________________________________________//
491
492 Float_t AliAnalysisKinkESDMC::GetSigmaToVertex(AliESDtrack* esdTrack) const {
493   // Calculates the number of sigma to the vertex.
494   
495   Float_t b[2];
496   Float_t bRes[2];
497   Float_t bCov[3];
498
499     esdTrack->GetImpactParameters(b,bCov);
500   
501   if (bCov[0]<=0 || bCov[2]<=0) {
502     //AliDebug(1, "Estimated b resolution lower or equal zero!");
503     bCov[0]=0; bCov[2]=0;
504   }
505   bRes[0] = TMath::Sqrt(bCov[0]);
506   bRes[1] = TMath::Sqrt(bCov[2]);
507   
508   if (bRes[0] == 0 || bRes[1] ==0) return -1;
509   
510   Float_t d = TMath::Sqrt(TMath::Power(b[0]/bRes[0],2) + TMath::Power(b[1]/bRes[1],2));
511   
512   if (TMath::Exp(-d * d / 2) < 1e-10) return 1000;
513   
514   d = TMath::ErfInverse(1 - TMath::Exp(-d * d / 2)) * TMath::Sqrt(2);
515   
516   return d;
517   }
518 //____________________________________________________________________//
519
520 const AliESDVertex* AliAnalysisKinkESDMC::GetEventVertex(const AliESDEvent* esd) const
521   // Get the vertex from the ESD and returns it if the vertex is valid
522   
523 {
524   // Get the vertex 
525   
526   const AliESDVertex* vertex = esd->GetPrimaryVertex();
527
528   if((vertex->GetStatus()==kTRUE)&&(vertex->GetNContributors()>2)) return vertex;
529   else
530   { 
531      vertex = esd->GetPrimaryVertexSPD();
532      if((vertex->GetStatus()==kTRUE)&&(vertex->GetNContributors()>2)) return vertex;
533      else
534      return 0;
535   }
536 }