]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/JetTasks/AliAnalyseLeadingTrackUE.cxx
Add histograms to check for special clusters with high energy and low number of cells
[u/mrichter/AliRoot.git] / PWG4 / JetTasks / AliAnalyseLeadingTrackUE.cxx
1 /*************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: A.Abrahantes, E.Lopez, S.Vallero                               *
5  * Version 1.0                                                            *
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 //#include <TBranch.h>
16 //#include <TCanvas.h>
17 //#include <TChain.h>
18 //#include <TFile.h>
19 //#include <TH1F.h>
20 //#include <TH1I.h>
21 //#include <TH2F.h>
22 #include <TList.h>
23 //#include <TLorentzVector.h>
24 #include <TMath.h>
25 #include <TObjArray.h>
26 #include <TObject.h>
27 //#include <TProfile.h>
28 //#include <TRandom.h>
29 //#include <TSystem.h>
30 //#include <TTree.h>
31 #include <TVector3.h>
32
33 #include "AliAnalyseLeadingTrackUE.h"
34 //#include "AliAnalysisTask.h"
35
36 //#include "AliAnalysisHelperJetTasks.h"
37 //#include "AliAnalysisManager.h"
38 #include "AliAODEvent.h"
39 //#include "AliAODHandler.h"
40 //#include "AliAODJet.h"
41 #include "AliAODMCParticle.h"
42 #include "AliAODTrack.h"
43 #include "AliESDEvent.h"
44 #include "AliESDtrack.h"
45 #include "AliESDtrackCuts.h"
46 //#include "AliGenPythiaEventHeader.h"
47 #include "AliInputEventHandler.h"
48 //#include "AliKFVertex.h"
49 //#include "AliLog.h"
50 #include "AliMCEvent.h"
51 //#include "AliMCEventHandler.h"
52 //#include "AliStack.h"
53 #include "AliVParticle.h"
54
55 ////////////////////////////////////////////////
56 //--------------------------------------------- 
57 // Class for transverse regions analysis
58 //---------------------------------------------
59 ////////////////////////////////////////////////
60
61
62 using namespace std;
63
64 ClassImp(AliAnalyseLeadingTrackUE)
65
66 //-------------------------------------------------------------------
67 AliAnalyseLeadingTrackUE::AliAnalyseLeadingTrackUE() :
68   TObject(),
69   fDebug(0),
70   fFilterBit(16),
71   fOnlyHadrons(kFALSE),
72   fTrackEtaCut(0.8),
73   fTrackPtMin(0),
74   fEsdTrackCuts(0x0), 
75   fEsdTrackCutsSPD(0x0), 
76   fEsdTrackCutsSDD(0x0) 
77 {
78   // constructor
79 }
80
81
82 //-------------------------------------------------------------------
83 AliAnalyseLeadingTrackUE & AliAnalyseLeadingTrackUE::operator = (const AliAnalyseLeadingTrackUE & /*source*/)
84 {
85   // assignment operator
86   return *this;
87 }
88
89
90 //-------------------------------------------------------------------
91 AliAnalyseLeadingTrackUE::~AliAnalyseLeadingTrackUE()
92 {
93
94   //clear memory
95   
96 }
97
98
99 //____________________________________________________________________
100 Bool_t AliAnalyseLeadingTrackUE::ApplyCuts(TObject* track)
101 {
102   
103   // select track according to set of cuts
104   if (! fEsdTrackCuts->IsSelected(track) )return kFALSE;
105   if (!fEsdTrackCutsSPD->IsSelected(track) && fEsdTrackCutsSDD->IsSelected(track)) return kFALSE;
106
107   return kTRUE;
108 }
109
110
111 //____________________________________________________________________
112 void AliAnalyseLeadingTrackUE::DefineESDCuts(Int_t /*filterbit*/){
113   
114    // Reproduces the cuts of the corresponding bit in the ESD->AOD filtering
115    // (see $ALICE_ROOT/ANALYSIS/macros/AddTaskESDFilter.C)
116
117    fEsdTrackCuts = AliESDtrackCuts::GetStandardITSTPCTrackCuts2010();
118    fEsdTrackCuts->SetClusterRequirementITS(AliESDtrackCuts::kSPD, AliESDtrackCuts::kNone);
119  
120    // Add SPD requirement 
121    fEsdTrackCutsSPD = new AliESDtrackCuts("SPD", "Require 1 cluster in SPD");
122    fEsdTrackCutsSPD->SetClusterRequirementITS(AliESDtrackCuts::kSPD,AliESDtrackCuts::kAny);
123  
124    // Add SDD requirement 
125    fEsdTrackCutsSDD = new AliESDtrackCuts("SDD", "Require 1 cluster in first layer SDD");
126    fEsdTrackCutsSDD->SetClusterRequirementITS(AliESDtrackCuts::kSDD,AliESDtrackCuts::kFirst);
127 }
128
129
130 //____________________________________________________________________
131 TObjArray*  AliAnalyseLeadingTrackUE::FindLeadingObjects(TObject *obj)
132 {
133
134   // Returns an array of charged particles (or jets) ordered according to their pT.
135
136   Int_t nTracks = NParticles(obj);
137
138
139   if( !nTracks ) return 0;
140  
141   // Define array of AliVParticle objects
142   TObjArray* tracks = new TObjArray(nTracks);
143
144   // Loop over tracks or jets
145   for (Int_t ipart=0; ipart<nTracks; ++ipart) {
146         AliVParticle* part = ParticleWithCuts( obj, ipart );
147         if (!part) continue;
148         // Accept leading-tracks in a limited pseudo-rapidity range     
149         if( TMath::Abs(part->Eta()) > fTrackEtaCut ) continue;
150         tracks->AddLast( part );
151         }
152   // Order tracks by pT 
153   QSortTracks( *tracks, 0, tracks->GetEntriesFast() );
154
155   nTracks = tracks->GetEntriesFast();
156   if( !nTracks ) return 0;
157
158   return tracks;
159   }
160
161
162 //-------------------------------------------------------------------
163 TObjArray* AliAnalyseLeadingTrackUE::GetAcceptedParticles(TObject* obj, TObject* arrayMC, Bool_t onlyprimaries, Int_t particleSpecies, Bool_t useEtaPtCuts)
164 {
165   // Returns an array of particles that pass the cuts, if arrayMC is given each reconstructed particle is replaced by its corresponding MC particles, depending on the parameter onlyprimaries only for primaries 
166   // particleSpecies: -1 all particles are returned
167   //                  0 (pions) 1 (kaons) 2 (protons) 3 (others) particles
168
169   Int_t nTracks = NParticles(obj);
170   TObjArray* tracks = new TObjArray;
171  
172   // Loop over tracks or jets
173   for (Int_t ipart=0; ipart<nTracks; ++ipart) {
174     AliVParticle* part = ParticleWithCuts( obj, ipart, onlyprimaries, particleSpecies );
175     if (!part) continue;
176     
177     if (useEtaPtCuts)
178       if (TMath::Abs(part->Eta()) > fTrackEtaCut || part->Pt() < fTrackPtMin)
179         continue;
180     
181     if (arrayMC && arrayMC->InheritsFrom("TClonesArray") && obj->InheritsFrom("AliAODEvent")) {
182       Int_t label = ((AliAODTrack*)part)->GetLabel();
183       // re-define part as the matched MC particle
184       part = (AliAODMCParticle*)ParticleWithCuts(arrayMC, TMath::Abs(label),onlyprimaries, particleSpecies);
185       if (!part)continue;
186     }
187     
188     tracks->Add(part);
189   }
190
191   return tracks;
192 }
193
194 //-------------------------------------------------------------------
195 TObjArray* AliAnalyseLeadingTrackUE::GetMinMaxRegion(TList *transv1, TList *transv2)
196 {
197   
198  // Returns two lists of particles, one for MIN and one for MAX region
199   Double_t sumpT1 = 0.;
200   Double_t sumpT2 = 0.;
201
202   Int_t particles1 = transv1->GetEntries();
203   Int_t particles2 = transv2->GetEntries();
204   
205   // Loop on transverse region 1
206   for (Int_t i=0; i<particles1; i++){
207         AliVParticle *part = (AliVParticle*)transv1->At(i);
208         sumpT1 +=  part->Pt();
209         }
210
211   // Loop on transverse region 2
212   for (Int_t i=0; i<particles2; i++){
213         AliVParticle *part = (AliVParticle*)transv2->At(i);
214         sumpT2 +=  part->Pt();
215         }
216
217   TObjArray *regionParticles = new TObjArray;
218   if ( sumpT2 >= sumpT1 ){
219         regionParticles->AddLast(transv1); // MIN
220         regionParticles->AddLast(transv2); // MAX 
221   }else {
222         regionParticles->AddLast(transv2); // MIN
223         regionParticles->AddLast(transv1); // MAX
224         }
225
226  return regionParticles;
227 }
228
229 //-------------------------------------------------------------------
230 Int_t  AliAnalyseLeadingTrackUE::NParticles(TObject* obj)
231 {
232  
233   //Returns the number of particles in AliAODMCParticle array  or AliAODTracks or AliESDTracks 
234
235   Int_t nTracks;
236   
237   if (obj->InheritsFrom("TClonesArray")){ // MC particles
238         TClonesArray *arrayMC = dynamic_cast<TClonesArray*>(obj);
239         nTracks = arrayMC->GetEntriesFast();
240   }else if (obj->InheritsFrom("TObjArray")){ // list of AliVParticle
241         TObjArray *array = dynamic_cast<TObjArray*>(obj);
242         nTracks = array->GetEntriesFast();
243   }else if (obj->InheritsFrom("AliAODEvent")){  // RECO AOD tracks
244         AliAODEvent *aodEvent = dynamic_cast<AliAODEvent*>(obj);
245         nTracks = aodEvent->GetNTracks();
246   }else if (obj->InheritsFrom("AliESDEvent")){  // RECO ESD tracks
247         AliESDEvent *esdEvent = dynamic_cast<AliESDEvent*>(obj);
248         nTracks = esdEvent->GetNumberOfTracks();
249   }else if (obj->InheritsFrom("AliMCEvent")){  // RECO ESD tracks
250         AliMCEvent *mcEvent = dynamic_cast<AliMCEvent*>(obj);
251         nTracks = mcEvent->GetNumberOfTracks();
252   }else {
253         if (fDebug > 1) AliFatal(" Analysis type not defined !!! ");
254         return 0;
255         }
256   
257   return nTracks;
258 }
259
260
261 //-------------------------------------------------------------------
262 AliVParticle*  AliAnalyseLeadingTrackUE::ParticleWithCuts(TObject* obj, Int_t ipart, Bool_t onlyprimaries, Int_t particleSpecies)
263 {
264   // Returns track or MC particle at position "ipart" if passes selection criteria
265   // particleSpecies: -1 all particles are returned
266   //                  0 (pions) 1 (kaons) 2 (protons) 3 (others) particles
267   AliVParticle *part=0;
268   
269   if (obj->InheritsFrom("TClonesArray")){ // AOD-MC PARTICLE
270         TClonesArray *arrayMC = dynamic_cast<TClonesArray*>(obj);
271         part = (AliVParticle*)arrayMC->At( ipart );
272         if (!part)return 0;
273         // eventually only primaries
274         if (onlyprimaries && !( ((AliAODMCParticle*)part)->IsPhysicalPrimary()) )return 0;
275         // eventually only hadrons
276         if (fOnlyHadrons){
277                 Int_t pdgCode = ((AliAODMCParticle*)part)->GetPdgCode();
278                 Bool_t isHadron = TMath::Abs(pdgCode)==211 ||  // Pion
279                                   TMath::Abs(pdgCode)==2212 || // Proton
280                                   TMath::Abs(pdgCode)==321;    // Kaon
281                 if (!isHadron) return 0;                                  
282                 }
283         if (particleSpecies != -1) {
284                 // find the primary mother
285                 AliVParticle* mother = part;
286                 while (!((AliAODMCParticle*)mother)->IsPhysicalPrimary())
287                 {
288                   if (((AliAODMCParticle*)mother)->GetMother() < 0)
289                   {
290                     mother = 0;
291                     break;
292                   }
293                     
294                   mother = (AliVParticle*) arrayMC->At(((AliAODMCParticle*)mother)->GetMother());
295                   if (!mother)
296                     break;
297                 }
298                 
299                 if (mother)
300                 {
301                   Int_t pdgCode = ((AliAODMCParticle*)mother)->GetPdgCode();
302                   if (particleSpecies == 0 && TMath::Abs(pdgCode)!=211)
303                           return 0;
304                   if (particleSpecies == 1 && TMath::Abs(pdgCode)!=321)
305                           return 0;
306                   if (particleSpecies == 2 && TMath::Abs(pdgCode)!=2212)
307                           return 0;
308                   if (particleSpecies == 3 && (TMath::Abs(pdgCode)==211 || TMath::Abs(pdgCode)==321 || TMath::Abs(pdgCode)==2212))
309                           return 0;
310                 }
311         }
312   
313   }else if (obj->InheritsFrom("TObjArray")){ // list of AliVParticle
314         TObjArray *array = dynamic_cast<TObjArray*>(obj);
315         part = (AliVParticle*)array->At( ipart );
316         if (!part)return 0;
317   }else if (obj->InheritsFrom("AliMCEvent")){ // MC PARTICLE
318         AliMCEvent* mcEvent =  dynamic_cast<AliMCEvent*>(obj);
319         part = mcEvent->GetTrack( ipart );
320         if (!part) return 0;
321         // eventually only primaries
322         if (onlyprimaries && !( mcEvent->IsPhysicalPrimary(ipart)) )return 0;
323         // eventually only hadrons
324         //TO-DO
325         /*if (fOnlyHadrons){
326                 Int_t pdgCode =  part->GetPdgCode();
327                 Bool_t isHadron = TMath::Abs(pdgCode)==211 ||  // Pion
328                                   TMath::Abs(pdgCode)==2212 || // Proton
329                                   TMath::Abs(pdgCode)==321;    // Kaon
330                 if (!isHadron) return 0;                                  
331                 }
332         */
333
334   }else if (obj->InheritsFrom("AliAODEvent")){ // RECO AOD TRACKS
335         AliAODEvent *aodEvent = dynamic_cast<AliAODEvent*>(obj);
336         part = aodEvent->GetTrack(ipart);
337         // track selection cuts
338         if ( !(((AliAODTrack*)part)->TestFilterBit(fFilterBit)) ) return 0; 
339         //if ( !(((AliAODTrack*)part)->TestFilterBit(fFilterBit)) && !(((AliAODTrack*)part)->TestFilterBit(32)) ) return 0; 
340         // only primary candidates
341         //if ( ((AliAODTrack*)part)->IsPrimaryCandidate() )return 0;
342         // eventually only hadrons
343         if (fOnlyHadrons){
344                 Bool_t isHadron = ((AliAODTrack*)part)->GetMostProbablePID()==AliAODTrack::kPion ||
345                                   ((AliAODTrack*)part)->GetMostProbablePID()==AliAODTrack::kKaon ||
346                                   ((AliAODTrack*)part)->GetMostProbablePID()==AliAODTrack::kProton;
347                 if (!isHadron) return 0;                                  
348                 }
349   
350   }else if (obj->InheritsFrom("AliESDEvent")){ // RECO ESD TRACKS
351         AliESDEvent *esdEvent = dynamic_cast<AliESDEvent*>(obj);
352         part = esdEvent->GetTrack(ipart);
353         if (!part)return 0;
354         // track selection cuts
355         
356         if (!( ApplyCuts(part)) )
357          return 0; 
358         
359         // only primary candidates (does not exist for ESD tracks??????)
360         //if ( ((AliAODTrack*)part)->IsPrimaryCandidate() )return 0;
361         
362         // eventually only hadrons
363         //TO-DO
364         /*if (fOnlyHadrons){
365                 Bool_t isHadron = ((AliAODTrack*)part)->GetMostProbablePID()==AliAODTrack::kPion ||
366                                   ((AliAODTrack*)part)->GetMostProbablePID()==AliAODTrack::kKaon ||
367                                   ((AliAODTrack*)part)->GetMostProbablePID()==AliAODTrack::kProton;
368                 if (!isHadron) return 0;                                  
369                 }
370         */      
371   }else {
372         if (fDebug > 1) AliFatal(" Analysis type not defined !!! ");
373         return 0;
374         }
375   
376   // only charged
377   if (!part->Charge())return 0;
378   
379   return part;
380 }
381
382
383 //-------------------------------------------------------------------
384 void  AliAnalyseLeadingTrackUE::QSortTracks(TObjArray &a, Int_t first, Int_t last)
385 {
386   // Sort array of TObjArray of tracks by Pt using a quicksort algorithm.
387   
388   static TObject *tmp;
389   static int i;           // "static" to save stack space
390   int j;
391   
392   while (last - first > 1) {
393     i = first;
394     j = last;
395     for (;;) {
396       while (++i < last && ((AliVParticle*)a[i])->Pt() > ((AliVParticle*)a[first])->Pt() )
397         ;
398       while (--j > first && ((AliVParticle*)a[j])->Pt() < ((AliVParticle*)a[first])->Pt() )
399         ;
400       if (i >= j)
401         break;
402       
403       tmp  = a[i];
404       a[i] = a[j];
405       a[j] = tmp;
406     }
407     if (j == first) {
408       ++first;
409       continue;
410     }
411     tmp = a[first];
412     a[first] = a[j];
413     a[j] = tmp;
414     if (j - first < last - (j + 1)) {
415       QSortTracks(a, first, j);
416       first = j + 1;   // QSortTracks(j + 1, last);
417     } else {
418       QSortTracks(a, j + 1, last);
419       last = j;        // QSortTracks(first, j);
420     }
421   }
422 }
423
424 //____________________________________________________________________
425 TObjArray*  AliAnalyseLeadingTrackUE::SortRegions(const AliVParticle* leading, TObject* obj, TObject* arrayMC, Bool_t onlyprimaries)
426 {
427
428   // Assign particles to towards, away or transverse regions.
429   // Returns a lists of particles for each region.
430
431   static const Double_t k60rad  = 60.*TMath::Pi()/180.;
432   static const Double_t k120rad = 120.*TMath::Pi()/180.;
433
434   // Define output lists of particles
435   TList *toward = new TList();
436   TList *away = new TList();
437   // Two transverse regions, for the moment those are not yet MIN and MAX!!! 
438   // MIN and MAX can be sorted in GetMinMaxRegion function
439   TList *transverse1 = new TList();
440   TList *transverse2 = new TList();
441   
442   TObjArray *regionParticles = new TObjArray;
443   regionParticles->SetOwner(kTRUE);
444   
445   regionParticles->AddLast(toward);
446   regionParticles->AddLast(away);
447   regionParticles->AddLast(transverse1);
448   regionParticles->AddLast(transverse2);
449   
450   if (!leading)
451     return regionParticles;
452  
453   // Switch to vector for leading particle
454   TVector3 leadVect(leading->Px(),leading->Py(),leading->Pz());
455   
456   Int_t nTracks = NParticles(obj);
457   if( !nTracks ) return 0;
458   // Loop over tracks 
459   for (Int_t ipart=0; ipart<nTracks; ++ipart) {
460         AliVParticle* part = ParticleWithCuts(obj, ipart);
461         if (!part)continue;
462         //Switch to vectors for particles 
463         TVector3 partVect(part->Px(), part->Py(), part->Pz());
464  
465         Int_t region = 0;
466         if( TMath::Abs(partVect.Eta()) > fTrackEtaCut ) continue;
467         // transverse regions
468         if (leadVect.DeltaPhi(partVect) < -k60rad && leadVect.DeltaPhi(partVect) > -k120rad )region = -1; //left
469         if (leadVect.DeltaPhi(partVect) > k60rad && leadVect.DeltaPhi(partVect) < k120rad ) region = 1;   //right
470
471         if (TMath::Abs(leadVect.DeltaPhi(partVect)) < k60rad ) region = 2;    //forward
472         if (TMath::Abs(leadVect.DeltaPhi(partVect)) > k120rad ) region = -2;  //backward
473         
474         // skip leading particle 
475         if (leading == part)
476           continue;
477         
478         if (!region)continue;
479         if (arrayMC && arrayMC->InheritsFrom("TClonesArray") && obj->InheritsFrom("AliAODEvent")){
480                 Int_t label = ((AliAODTrack*)part)->GetLabel();
481                 // re-define part as the matched MC particle
482                 part = (AliAODMCParticle*)ParticleWithCuts(arrayMC, TMath::Abs(label),onlyprimaries);
483                 if (!part)continue;
484                 // skip leading particle 
485                 if (leading == part)
486                   continue;
487                 }
488         if (arrayMC && arrayMC->InheritsFrom("AliMCEvent") && obj->InheritsFrom("AliESDEvent")){
489                 Int_t label = ((AliESDtrack*)part)->GetLabel();
490                 // look for the matched MC particle (but do not re-define part)
491                 if (!ParticleWithCuts(arrayMC, TMath::Abs(label),onlyprimaries)) continue;
492                 }
493
494         if ( region == 1 ) transverse1->Add(part);
495         if ( region == -1 ) transverse2->Add(part);
496         if ( region == 2 ) toward->Add(part);
497         if ( region == -2 ) away->Add(part);
498
499         }//end loop on tracks
500   
501   return regionParticles;
502   
503 }
504
505
506 //____________________________________________________________________
507 Bool_t  AliAnalyseLeadingTrackUE::TriggerSelection(const TObject* obj)
508 {
509
510   //Use AliPhysicsSelection to select good events
511   if( !obj->InheritsFrom("AliAODInputHandler") ) { // not needed for AOD input 
512         if (!(((AliInputEventHandler*)obj)->IsEventSelected()&AliVEvent::kMB))return kFALSE;
513         }                                
514
515   // TODO for AOD input trigger bit needs to be checked too (is stored in the AOD)
516
517   return kTRUE;
518
519 }
520
521 //____________________________________________________________________
522 Bool_t  AliAnalyseLeadingTrackUE::VertexSelection(const TObject* obj, Int_t ntracks, Double_t zed)
523 {
524
525   //Require 1 vertex (no TPC stand-alone) with a minimum number of tracks and z-coordinate in a limited range
526   
527   if (obj->InheritsFrom("AliAODEvent")){ 
528         Int_t nVertex = ((AliAODEvent*)obj)->GetNumberOfVertices();
529         if( nVertex > 0 ) { 
530                 AliAODVertex* vertex = (AliAODVertex*)((AliAODEvent*)obj)->GetPrimaryVertex();
531                 Int_t nTracksPrim = vertex->GetNContributors();
532                 Double_t zVertex = vertex->GetZ();
533                 if (fDebug > 1)AliInfo(Form(" Vertex in = %f with %d particles by  %s data ...",zVertex,nTracksPrim,vertex->GetName()));
534                 // Reject TPC only vertex
535                 TString name(vertex->GetName());
536                 if (name.CompareTo("PrimaryVertex") && name.CompareTo("SPDVertex"))return kFALSE;
537
538                 // Select a quality vertex by number of tracks?
539                 if( nTracksPrim < ntracks || TMath::Abs(zVertex) > zed ) {
540                         if (fDebug > 1) AliInfo(" Primary-vertex Selection: event REJECTED ...");
541                         return kFALSE;
542                         }
543                 // TODO remove vertexer Z events with dispersion > 0.02: Doesn't work for AOD at present
544                 //if (strcmp(vertex->GetTitle(), "AliVertexerZ") == 0 && vertex->GetDispersion() > 0.02)
545                 //  return kFALSE;
546                 if (fDebug > 1) AliInfo(" Primary-vertex Selection: event ACCEPTED...");
547                 } else {
548                         if (fDebug > 1) AliInfo(" Primary-vertex Selection: event REJECTED ...");
549                         return kFALSE;
550                         }
551         }
552
553   if (obj->InheritsFrom("AliMCEvent"))
554   { 
555     if (TMath::Abs(((AliMCEvent*) obj)->GetPrimaryVertex()->GetZ()) > zed)
556     {
557       if (fDebug > 1) AliInfo(" Primary-vertex Selection: event (based on MC) REJECTED ...");
558       return kFALSE;
559     }
560   }
561
562   // ESD case for DCA studies
563   if (obj->InheritsFrom("AliESDEvent")){
564        AliESDVertex* vertex = (AliESDVertex*)((AliESDEvent*)obj)->GetPrimaryVertex();
565        if ( vertex){
566                Int_t nTracksPrim = vertex->GetNContributors();
567                Double_t zVertex = vertex->GetZ();
568                if (fDebug > 1)AliInfo(Form(" Vertex in = %f with %d particles by  %s data ...",zVertex,nTracksPrim,vertex->GetName()));
569                // Reject SPD or TPC only vertex
570                TString name(vertex->GetName());
571                if (name.CompareTo("PrimaryVertex") && name.CompareTo("SPDVertex"))return kFALSE;
572
573                // Select a quality vertex by number of tracks?
574                if( nTracksPrim < ntracks || TMath::Abs(zVertex) > zed ) {
575                        if (fDebug > 1) AliInfo(" Primary-vertex Selection: event REJECTED ...");
576                        return kFALSE;
577                        }
578                // TODO remove vertexer Z events with dispersion > 0.02: Doesn't work for AOD at present
579                 //if (strcmp(vertex->GetTitle(), "AliVertexerZ") == 0 && vertex->GetDispersion() > 0.02)
580                 //  return kFALSE;
581                if (fDebug > 1) AliInfo(" Primary-vertex Selection: event ACCEPTED...");
582                } else {
583                        if (fDebug > 1) AliInfo(" Primary-vertex Selection: event REJECTED ...");
584                        return kFALSE;
585                        }
586        }
587         
588   return kTRUE;
589 }