]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/KINK/AliAnalysisKinkESDMC.cxx
Newlines at EOF
[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(fAngMomK);
187    fListOfHistos->Add(fAngMomPi);
188    fListOfHistos->Add(fRpr);
189    fListOfHistos->Add(fZpr);
190    fListOfHistos->Add(fAngMomKC);
191    fListOfHistos->Add(fZvXv);
192    fListOfHistos->Add(fZvYv);
193    fListOfHistos->Add(fXvYv);
194    fListOfHistos->Add(fPtPrKink);
195 }
196
197 //________________________________________________________________________
198 void AliAnalysisKinkESDMC::Exec(Option_t *) 
199 {
200   // Main loop
201   // Called for each event
202
203   // Process MC truth, therefore we receive the AliAnalysisManager and ask it for the AliMCEventHandler
204   // This handler can return the current MC event
205  
206   AliMCEventHandler* eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
207   if (!eventHandler) {
208     Printf("ERROR: Could not retrieve MC event handler");
209     return;
210   }
211
212   AliMCEvent* mcEvent = eventHandler->MCEvent();
213   if (!mcEvent) {
214      Printf("ERROR: Could not retrieve MC event");
215      return;
216   }
217
218   Printf("MC particles: %d", mcEvent->GetNumberOfTracks());
219
220   if (!fESD) {
221     Printf("ERROR: fESD not available");
222     return;
223   }
224   AliStack* stack=mcEvent->Stack();
225   
226 //primary tracks  in MC
227          Int_t  nPrim = stack->GetNprimary();
228 //
229   const AliESDVertex *vertex=GetEventVertex(fESD);
230     if(!vertex) return;
231 //
232   Double_t vpos[3];
233   vertex->GetXYZ(vpos);
234     fZpr->Fill(vpos[2]);         
235  fZvXv->Fill( vpos[1], vpos[2]);  
236  fZvYv->Fill( vpos[0], vpos[2]);  
237
238   Double_t vtrack[3], ptrack[3];
239   
240      
241
242    Int_t nGoodTracks =  fESD->GetNumberOfTracks();
243     fESDMult->Fill(nGoodTracks);
244 //
245 // track loop
246   for (Int_t iTracks = 0; iTracks < fESD->GetNumberOfTracks(); iTracks++) {
247     AliESDtrack* track = fESD->GetTrack(iTracks);
248     if (!track) {
249       Printf("ERROR: Could not receive track %d", iTracks);
250       continue;
251     }
252     
253
254     fHistPt->Fill(track->Pt());
255
256       //   pt cut at 300 MeV
257         if ( (track->Pt())<.300)continue;
258
259     UInt_t status=track->GetStatus();
260
261   //  if((status&AliESDtrack::kITSrefit)==0) continue;
262     if((status&AliESDtrack::kTPCrefit)==0) continue;
263       if((track->GetTPCchi2()/track->GetTPCclusters(0))>3.5) continue;
264
265       Double_t extCovPos[15];
266       track->GetExternalCovariance(extCovPos);    
267       if(extCovPos[0]>2) continue;
268       if(extCovPos[2]>2) continue;    
269       if(extCovPos[5]>0.5) continue;  
270       if(extCovPos[9]>0.5) continue;
271       if(extCovPos[14]>2) continue;
272
273
274     track->GetXYZ(vtrack);
275  fXvYv->Fill(vtrack[0],vtrack[1]);  
276
277 // track momentum
278      track->GetPxPyPz(ptrack);
279     
280     TVector3 trackMom(ptrack[0],ptrack[1],ptrack[2]);
281     
282     Double_t trackEta=trackMom.Eta();
283     
284     
285     
286     Float_t nSigmaToVertex = GetSigmaToVertex(track);  
287       
288     Float_t bpos[2];
289     Float_t bCovpos[3];
290     track->GetImpactParameters(bpos,bCovpos);
291     
292     if (bCovpos[0]<=0 || bCovpos[2]<=0) {
293      Printf("Estimated b resolution lower or equal zero!");
294      bCovpos[0]=0; bCovpos[2]=0;
295     }
296
297     Float_t dcaToVertexXYpos = bpos[0];
298     Float_t dcaToVertexZpos = bpos[1];
299     
300     fRpr->Fill(dcaToVertexZpos);
301
302     if(nSigmaToVertex>=4) continue;
303     if((dcaToVertexXYpos>3.0)||(dcaToVertexZpos>3.0)) continue;
304     
305     
306     Int_t label = track->GetLabel();
307     label = TMath::Abs(label);
308     TParticle * part = stack->Particle(label);
309     if (!part) continue;
310
311  //  cut on eta 
312         if(  (TMath::Abs(trackEta )) > 1.1 ) continue;
313     fHistPtESD->Fill(track->Pt());
314
315    // Add Kink analysis
316    
317                 Int_t indexKinkPos=track->GetKinkIndex(0);
318 //  loop on kinks
319                 if(indexKinkPos<0){
320                fPtKink->Fill(track->Pt()); /// pt from track
321
322         // select kink class    
323
324           AliESDkink *kink=fESD->GetKink(TMath::Abs(indexKinkPos)-1);
325 //
326         
327           Int_t eSDfLabel1=kink->GetLabel(0);
328           TParticle *particle1= stack->Particle(TMath::Abs(eSDfLabel1));
329           Int_t code1= particle1->GetPdgCode();
330           
331           Int_t eSDfLabeld=kink->GetLabel(1);
332           TParticle *particled= stack->Particle(TMath::Abs(eSDfLabeld));
333           Int_t dcode1= particled->GetPdgCode();
334
335            const TVector3 motherMfromKink(kink->GetMotherP());
336            const TVector3 daughterMKink(kink->GetDaughterP());
337            Float_t qT=kink->GetQt();
338
339            fHistEta->Fill(trackEta) ;  //   Eta distr
340            fHistQtAll->Fill(qT) ;  //  Qt   distr
341                   
342            fptKink->Fill(motherMfromKink.Pt()); /// pt from kink
343
344 //       fake kinks, small Qt and small kink angle
345            if ( qT<0.012) fcodeH->Fill( TMath::Abs(code1), TMath::Abs(dcode1) );
346            if(qT<0.012) continue;  // remove fake kinks
347
348            Float_t kinkAngle=TMath::RadToDeg()*kink->GetAngle(2);
349          
350 //remove the double taracks 
351             
352            if( (kinkAngle<1.)  ) continue;
353 //
354            if( (qT>0.05)&& ( qT<0.25)  )   fHistQt2->Fill(qT);  // candidate kaon kinks
355
356            if (TMath::Abs(code1)==321)     fHistPtKPDG->Fill(track->Pt());  // ALL KAONS (pdg) inside kink sample
357
358 //          maximum decay angle at a given mother momentum
359            Double_t maxDecAngKmu=f1->Eval(motherMfromKink.Mag(),0.,0.,0.);
360            Double_t maxDecAngpimu=f2->Eval(motherMfromKink.Mag(),0.,0.,0.);
361 //  two dimensional plot 
362                 if(TMath::Abs(code1)==321) fAngMomK->Fill(motherMfromKink.Mag(), kinkAngle); 
363                 if(TMath::Abs(code1)==211) fAngMomPi->Fill(motherMfromKink.Mag(), kinkAngle); 
364 //
365 // invariant mass of mother track decaying to mu
366          Float_t energyDaughterMu=TMath::Sqrt(daughterMKink.Mag()*daughterMKink.Mag()+0.105658*0.105658);
367          Float_t p1XM= motherMfromKink.Px();
368          Float_t p1YM= motherMfromKink.Py();
369          Float_t p1ZM= motherMfromKink.Pz();
370          Float_t p2XM= daughterMKink.Px();
371          Float_t p2YM= daughterMKink.Py();
372          Float_t p2ZM= daughterMKink.Pz();
373          Float_t p3Daughter=TMath::Sqrt(((p1XM-p2XM)*(p1XM-p2XM))+((p1YM-p2YM)*(p1YM-p2YM))+((p1ZM-p2ZM)*(p1ZM-p2ZM)));
374          Double_t invariantMassKmu= TMath::Sqrt((energyDaughterMu+p3Daughter)*(energyDaughterMu+p3Daughter)-motherMfromKink.Mag()*motherMfromKink.Mag());
375   
376          fM1kaon->Fill(invariantMassKmu);
377
378 //  kaon selection from kinks
379 if((kinkAngle>maxDecAngpimu)&&(qT>0.05)&&(qT<0.25)&&((kink->GetR()>110.)&&(kink->GetR()<230.))&&(TMath::Abs(trackEta)<1.1)&&(invariantMassKmu<0.6)) {
380                 fHistEtaK->Fill(trackEta);
381
382                 if(TMath::Abs(code1)==321) fAngMomKC->Fill(motherMfromKink.Mag(), kinkAngle);  // real kaons
383
384                 fHistQt1  ->Fill(qT) ;  //  Qt   distr
385                 if ( (kinkAngle>maxDecAngKmu) && ( motherMfromKink.Mag()> 1.1 ) ) continue; // maximum angle selection
386
387                                     fHistPtKaon->Fill(track->Pt());   //all PID kink-kaon
388
389 // background inside the identified kaons, e.g KK  ,  Kp
390          if((TMath::Abs(code1)==321)&&(( TMath::Abs(dcode1)) >=(TMath::Abs(code1)) )) fdcodeH->Fill( TMath::Abs(code1), TMath::Abs(dcode1));
391
392          if((TMath::Abs(code1)==321) ) { 
393
394          if ( label< nPrim ) fPtPrKink->Fill(track->Pt());
395              fKinkKaon->Fill(track->Pt());        }
396          else {
397              fKinkKaonBg->Fill(track->Pt());     
398           }   // primary and secondary 
399
400         }  //  kink selection 
401                   
402
403         }  //End Kink Information    
404   
405
406   } //track loop 
407
408   // loop over mc particles
409
410   // variable to count tracks
411   Int_t nMCTracks = 0;
412
413 for (Int_t iMc = 0; iMc < nPrim; ++iMc)
414   {
415    // AliDebug(AliLog::kDebug+1, Form("MC Loop: Processing particle %d.", iMc));
416
417     TParticle* particle = stack->Particle(iMc);
418
419     if (!particle)
420     {
421     //  AliDebug(AliLog::kError, Form("UNEXPECTED: particle with label %d not found in stack (mc loop).", iMc));
422       continue;
423     }
424
425  //   
426     Float_t eta = particle->Eta();
427
428     if (TMath::Abs(eta) > 1.1)
429       continue;
430
431             Double_t ptK = particle->Pt();
432
433    if( ptK <0.300) continue;
434
435       Float_t code = particle->GetPdgCode();
436       if (code==321 || code==-321){
437             
438               fptKMC->Fill(ptK);
439     
440
441             Int_t firstD=particle->GetFirstDaughter();
442             Int_t lastD=particle->GetLastDaughter();
443 //loop on secondaries
444              for (Int_t k=firstD;k<lastD;k++) {
445              TParticle*    daughter1=stack->Particle(k+1);
446              Float_t dcode = daughter1->GetPdgCode();
447
448     
449
450                  frad->Fill(daughter1->R());
451
452                  if (((daughter1->R())>110)&&((daughter1->R())<230) ){
453         if (( ( code==321 )&& ( dcode ==-13  ))|| (( code == -321 )&& ( dcode == 13))) fgenPtEtR->Fill( ptK );//to muon
454         if (( ( code==321 )&& ( dcode ==211  ))|| (( code == -321 )&& ( dcode ==-211))) fgenPtEtR->Fill( ptK );//to pion
455      //   if (( ( code==321 )&& ( dcode ==-11  ))|| (( code == -321 )&& ( dcode ==11))) fgenPtEtR->Fill( ptK );// to electr
456
457        }
458 //  next only  to mu decay
459                  if (((code==321)&&(dcode==-13))||((code==-321)&&(dcode==13))){
460
461                  
462                     if (((daughter1->R())>110)&&((daughter1->R())<230)){
463                     fgenpt->Fill(ptK);
464                     }
465                     }
466        }
467
468       }
469
470
471       
472     nMCTracks++;
473   }// end of mc particle
474
475 //  printf("hello mc");
476   fMultiplMC->Fill(nMCTracks);
477
478   PostData(1, fListOfHistos);
479
480 }      
481
482 //________________________________________________________________________
483 void AliAnalysisKinkESDMC::Terminate(Option_t *) 
484 {
485   // Draw result to the screen 
486   // Called once at the end of the query
487
488 }
489 //____________________________________________________________________//
490
491 Float_t AliAnalysisKinkESDMC::GetSigmaToVertex(AliESDtrack* esdTrack) const {
492   // Calculates the number of sigma to the vertex.
493   
494   Float_t b[2];
495   Float_t bRes[2];
496   Float_t bCov[3];
497
498     esdTrack->GetImpactParameters(b,bCov);
499   
500   if (bCov[0]<=0 || bCov[2]<=0) {
501     //AliDebug(1, "Estimated b resolution lower or equal zero!");
502     bCov[0]=0; bCov[2]=0;
503   }
504   bRes[0] = TMath::Sqrt(bCov[0]);
505   bRes[1] = TMath::Sqrt(bCov[2]);
506   
507   if (bRes[0] == 0 || bRes[1] ==0) return -1;
508   
509   Float_t d = TMath::Sqrt(TMath::Power(b[0]/bRes[0],2) + TMath::Power(b[1]/bRes[1],2));
510   
511   if (TMath::Exp(-d * d / 2) < 1e-10) return 1000;
512   
513   d = TMath::ErfInverse(1 - TMath::Exp(-d * d / 2)) * TMath::Sqrt(2);
514   
515   return d;
516   }
517 //____________________________________________________________________//
518
519 const AliESDVertex* AliAnalysisKinkESDMC::GetEventVertex(const AliESDEvent* esd) const
520   // Get the vertex from the ESD and returns it if the vertex is valid
521   
522 {
523   // Get the vertex 
524   
525   const AliESDVertex* vertex = esd->GetPrimaryVertex();
526
527   if((vertex->GetStatus()==kTRUE)&&(vertex->GetNContributors()>2)) return vertex;
528   else
529   { 
530      vertex = esd->GetPrimaryVertexSPD();
531      if((vertex->GetStatus()==kTRUE)&&(vertex->GetNContributors()>2)) return vertex;
532      else
533      return 0;
534   }
535 }