]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/JetTasks/AliAnalyseLeadingTrackUE.cxx
run on AOD branches
[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                 else
312                 {
313                   // if mother not found, accept particle only in case of particleSpecies == 3. To include it in all or no sample is no solution
314                   Printf("WARNING: No mother found for particle %d", part->GetLabel());
315                   if (particleSpecies != 3)
316                     return 0;
317                 }
318         }
319   
320   }else if (obj->InheritsFrom("TObjArray")){ // list of AliVParticle
321         TObjArray *array = dynamic_cast<TObjArray*>(obj);
322         part = (AliVParticle*)array->At( ipart );
323         if (!part)return 0;
324   }else if (obj->InheritsFrom("AliMCEvent")){ // MC PARTICLE
325         AliMCEvent* mcEvent =  dynamic_cast<AliMCEvent*>(obj);
326         part = mcEvent->GetTrack( ipart );
327         if (!part) return 0;
328         // eventually only primaries
329         if (onlyprimaries && !( mcEvent->IsPhysicalPrimary(ipart)) )return 0;
330         // eventually only hadrons
331         //TO-DO
332         /*if (fOnlyHadrons){
333                 Int_t pdgCode =  part->GetPdgCode();
334                 Bool_t isHadron = TMath::Abs(pdgCode)==211 ||  // Pion
335                                   TMath::Abs(pdgCode)==2212 || // Proton
336                                   TMath::Abs(pdgCode)==321;    // Kaon
337                 if (!isHadron) return 0;                                  
338                 }
339         */
340
341   }else if (obj->InheritsFrom("AliAODEvent")){ // RECO AOD TRACKS
342         AliAODEvent *aodEvent = dynamic_cast<AliAODEvent*>(obj);
343         part = aodEvent->GetTrack(ipart);
344         // track selection cuts
345         if ( !(((AliAODTrack*)part)->TestFilterBit(fFilterBit)) ) return 0; 
346         //if ( !(((AliAODTrack*)part)->TestFilterBit(fFilterBit)) && !(((AliAODTrack*)part)->TestFilterBit(32)) ) return 0; 
347         // only primary candidates
348         //if ( ((AliAODTrack*)part)->IsPrimaryCandidate() )return 0;
349         // eventually only hadrons
350         if (fOnlyHadrons){
351                 Bool_t isHadron = ((AliAODTrack*)part)->GetMostProbablePID()==AliAODTrack::kPion ||
352                                   ((AliAODTrack*)part)->GetMostProbablePID()==AliAODTrack::kKaon ||
353                                   ((AliAODTrack*)part)->GetMostProbablePID()==AliAODTrack::kProton;
354                 if (!isHadron) return 0;                                  
355                 }
356   
357   }else if (obj->InheritsFrom("AliESDEvent")){ // RECO ESD TRACKS
358         AliESDEvent *esdEvent = dynamic_cast<AliESDEvent*>(obj);
359         part = esdEvent->GetTrack(ipart);
360         if (!part)return 0;
361         // track selection cuts
362         
363         if (!( ApplyCuts(part)) )
364          return 0; 
365         
366         // only primary candidates (does not exist for ESD tracks??????)
367         //if ( ((AliAODTrack*)part)->IsPrimaryCandidate() )return 0;
368         
369         // eventually only hadrons
370         //TO-DO
371         /*if (fOnlyHadrons){
372                 Bool_t isHadron = ((AliAODTrack*)part)->GetMostProbablePID()==AliAODTrack::kPion ||
373                                   ((AliAODTrack*)part)->GetMostProbablePID()==AliAODTrack::kKaon ||
374                                   ((AliAODTrack*)part)->GetMostProbablePID()==AliAODTrack::kProton;
375                 if (!isHadron) return 0;                                  
376                 }
377         */      
378   }else {
379         if (fDebug > 1) AliFatal(" Analysis type not defined !!! ");
380         return 0;
381         }
382   
383   // only charged
384   if (!part->Charge())return 0;
385   
386   return part;
387 }
388
389
390 //-------------------------------------------------------------------
391 void  AliAnalyseLeadingTrackUE::QSortTracks(TObjArray &a, Int_t first, Int_t last)
392 {
393   // Sort array of TObjArray of tracks by Pt using a quicksort algorithm.
394   
395   static TObject *tmp;
396   static int i;           // "static" to save stack space
397   int j;
398   
399   while (last - first > 1) {
400     i = first;
401     j = last;
402     for (;;) {
403       while (++i < last && ((AliVParticle*)a[i])->Pt() > ((AliVParticle*)a[first])->Pt() )
404         ;
405       while (--j > first && ((AliVParticle*)a[j])->Pt() < ((AliVParticle*)a[first])->Pt() )
406         ;
407       if (i >= j)
408         break;
409       
410       tmp  = a[i];
411       a[i] = a[j];
412       a[j] = tmp;
413     }
414     if (j == first) {
415       ++first;
416       continue;
417     }
418     tmp = a[first];
419     a[first] = a[j];
420     a[j] = tmp;
421     if (j - first < last - (j + 1)) {
422       QSortTracks(a, first, j);
423       first = j + 1;   // QSortTracks(j + 1, last);
424     } else {
425       QSortTracks(a, j + 1, last);
426       last = j;        // QSortTracks(first, j);
427     }
428   }
429 }
430
431 //____________________________________________________________________
432 TObjArray*  AliAnalyseLeadingTrackUE::SortRegions(const AliVParticle* leading, TObject* obj, TObject* arrayMC, Bool_t onlyprimaries)
433 {
434
435   // Assign particles to towards, away or transverse regions.
436   // Returns a lists of particles for each region.
437
438   static const Double_t k60rad  = 60.*TMath::Pi()/180.;
439   static const Double_t k120rad = 120.*TMath::Pi()/180.;
440
441   // Define output lists of particles
442   TList *toward = new TList();
443   TList *away = new TList();
444   // Two transverse regions, for the moment those are not yet MIN and MAX!!! 
445   // MIN and MAX can be sorted in GetMinMaxRegion function
446   TList *transverse1 = new TList();
447   TList *transverse2 = new TList();
448   
449   TObjArray *regionParticles = new TObjArray;
450   regionParticles->SetOwner(kTRUE);
451   
452   regionParticles->AddLast(toward);
453   regionParticles->AddLast(away);
454   regionParticles->AddLast(transverse1);
455   regionParticles->AddLast(transverse2);
456   
457   if (!leading)
458     return regionParticles;
459  
460   // Switch to vector for leading particle
461   TVector3 leadVect(leading->Px(),leading->Py(),leading->Pz());
462   
463   Int_t nTracks = NParticles(obj);
464   if( !nTracks ) return 0;
465   // Loop over tracks 
466   for (Int_t ipart=0; ipart<nTracks; ++ipart) {
467         AliVParticle* part = ParticleWithCuts(obj, ipart);
468         if (!part)continue;
469         //Switch to vectors for particles 
470         TVector3 partVect(part->Px(), part->Py(), part->Pz());
471  
472         Int_t region = 0;
473         if( TMath::Abs(partVect.Eta()) > fTrackEtaCut ) continue;
474         // transverse regions
475         if (leadVect.DeltaPhi(partVect) < -k60rad && leadVect.DeltaPhi(partVect) > -k120rad )region = -1; //left
476         if (leadVect.DeltaPhi(partVect) > k60rad && leadVect.DeltaPhi(partVect) < k120rad ) region = 1;   //right
477
478         if (TMath::Abs(leadVect.DeltaPhi(partVect)) < k60rad ) region = 2;    //forward
479         if (TMath::Abs(leadVect.DeltaPhi(partVect)) > k120rad ) region = -2;  //backward
480         
481         // skip leading particle 
482         if (leading == part)
483           continue;
484         
485         if (!region)continue;
486         if (arrayMC && arrayMC->InheritsFrom("TClonesArray") && obj->InheritsFrom("AliAODEvent")){
487                 Int_t label = ((AliAODTrack*)part)->GetLabel();
488                 // re-define part as the matched MC particle
489                 part = (AliAODMCParticle*)ParticleWithCuts(arrayMC, TMath::Abs(label),onlyprimaries);
490                 if (!part)continue;
491                 // skip leading particle 
492                 if (leading == part)
493                   continue;
494                 }
495         if (arrayMC && arrayMC->InheritsFrom("AliMCEvent") && obj->InheritsFrom("AliESDEvent")){
496                 Int_t label = ((AliESDtrack*)part)->GetLabel();
497                 // look for the matched MC particle (but do not re-define part)
498                 if (!ParticleWithCuts(arrayMC, TMath::Abs(label),onlyprimaries)) continue;
499                 }
500
501         if ( region == 1 ) transverse1->Add(part);
502         if ( region == -1 ) transverse2->Add(part);
503         if ( region == 2 ) toward->Add(part);
504         if ( region == -2 ) away->Add(part);
505
506         }//end loop on tracks
507   
508   return regionParticles;
509   
510 }
511
512
513 //____________________________________________________________________
514 Bool_t  AliAnalyseLeadingTrackUE::TriggerSelection(const TObject* obj)
515 {
516   if (!obj) // MC
517     return kFALSE;
518
519   //Use AliPhysicsSelection to select good events
520   if( !obj->InheritsFrom("AliAODInputHandler") ) { // not needed for AOD input 
521         if (!(((AliInputEventHandler*)obj)->IsEventSelected()&AliVEvent::kMB))return kFALSE;
522         }                                
523
524   // TODO for AOD input trigger bit needs to be checked too (is stored in the AOD)
525
526   return kTRUE;
527 }
528
529 //____________________________________________________________________
530 Bool_t  AliAnalyseLeadingTrackUE::VertexSelection(const TObject* obj, Int_t ntracks, Double_t zed)
531 {
532
533   //Require 1 vertex (no TPC stand-alone) with a minimum number of tracks and z-coordinate in a limited range
534   
535   if (obj->InheritsFrom("AliAODEvent")){ 
536         Int_t nVertex = ((AliAODEvent*)obj)->GetNumberOfVertices();
537         if( nVertex > 0 ) { 
538                 AliAODVertex* vertex = (AliAODVertex*)((AliAODEvent*)obj)->GetPrimaryVertex();
539                 Int_t nTracksPrim = vertex->GetNContributors();
540                 Double_t zVertex = vertex->GetZ();
541                 if (fDebug > 1)AliInfo(Form(" Vertex in = %f with %d particles by  %s data ...",zVertex,nTracksPrim,vertex->GetName()));
542                 // Reject TPC only vertex
543                 TString name(vertex->GetName());
544                 if (name.CompareTo("PrimaryVertex") && name.CompareTo("SPDVertex"))return kFALSE;
545
546                 // Select a quality vertex by number of tracks?
547                 if( nTracksPrim < ntracks || TMath::Abs(zVertex) > zed ) {
548                         if (fDebug > 1) AliInfo(" Primary-vertex Selection: event REJECTED ...");
549                         return kFALSE;
550                         }
551                 // TODO remove vertexer Z events with dispersion > 0.02: Doesn't work for AOD at present
552                 //if (strcmp(vertex->GetTitle(), "AliVertexerZ") == 0 && vertex->GetDispersion() > 0.02)
553                 //  return kFALSE;
554                 if (fDebug > 1) AliInfo(" Primary-vertex Selection: event ACCEPTED...");
555                 } else {
556                         if (fDebug > 1) AliInfo(" Primary-vertex Selection: event REJECTED ...");
557                         return kFALSE;
558                         }
559         }
560
561   if (obj->InheritsFrom("AliMCEvent"))
562   { 
563     if (TMath::Abs(((AliMCEvent*) obj)->GetPrimaryVertex()->GetZ()) > zed)
564     {
565       if (fDebug > 1) AliInfo(" Primary-vertex Selection: event (based on MC) REJECTED ...");
566       return kFALSE;
567     }
568   }
569
570   // ESD case for DCA studies
571   if (obj->InheritsFrom("AliESDEvent")){
572        AliESDVertex* vertex = (AliESDVertex*)((AliESDEvent*)obj)->GetPrimaryVertex();
573        if ( vertex){
574                Int_t nTracksPrim = vertex->GetNContributors();
575                Double_t zVertex = vertex->GetZ();
576                if (fDebug > 1)AliInfo(Form(" Vertex in = %f with %d particles by  %s data ...",zVertex,nTracksPrim,vertex->GetName()));
577                // Reject SPD or TPC only vertex
578                TString name(vertex->GetName());
579                if (name.CompareTo("PrimaryVertex") && name.CompareTo("SPDVertex"))return kFALSE;
580
581                // Select a quality vertex by number of tracks?
582                if( nTracksPrim < ntracks || TMath::Abs(zVertex) > zed ) {
583                        if (fDebug > 1) AliInfo(" Primary-vertex Selection: event REJECTED ...");
584                        return kFALSE;
585                        }
586                // TODO remove vertexer Z events with dispersion > 0.02: Doesn't work for AOD at present
587                 //if (strcmp(vertex->GetTitle(), "AliVertexerZ") == 0 && vertex->GetDispersion() > 0.02)
588                 //  return kFALSE;
589                if (fDebug > 1) AliInfo(" Primary-vertex Selection: event ACCEPTED...");
590                } else {
591                        if (fDebug > 1) AliInfo(" Primary-vertex Selection: event REJECTED ...");
592                        return kFALSE;
593                        }
594        }
595         
596   return kTRUE;
597 }