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