]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliAODEvent.cxx
Pile-up vertices added to AODEvent. (F. Prino)
[u/mrichter/AliRoot.git] / STEER / AliAODEvent.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2007, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
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 /* $Id$ */
17
18 //-------------------------------------------------------------------------
19 //     AOD base class
20 //     Author: Markus Oldenburg, CERN
21 //-------------------------------------------------------------------------
22
23 #include <TROOT.h>
24 #include <TTree.h>
25 #include <TFolder.h>
26 #include <TFriendElement.h>
27 #include <TProcessID.h>
28 #include <TCollection.h>
29 #include "Riostream.h"
30 #include "AliAODEvent.h"
31 #include "AliAODHeader.h"
32 #include "AliAODTrack.h"
33 #include "AliAODDimuon.h"
34
35 ClassImp(AliAODEvent)
36
37 // definition of std AOD member names
38   const char* AliAODEvent::fAODListName[kAODListN] = {"header",
39                                                       "tracks",
40                                                       "vertices",
41                                                       "v0s",
42                                                       "cascades",
43                                                       "tracklets",
44                                                       "jets",
45                                                       "emcalCells",
46                                                       "phosCells",
47                                                       "caloClusters",
48                                                       "fmdClusters",
49                                                       "pmdClusters",
50                                                       "dimuons"
51                                                       
52 };
53 //______________________________________________________________________________
54 AliAODEvent::AliAODEvent() :
55   AliVEvent(),
56   fAODObjects(new TList()),
57   fAODFolder(0),
58   fConnected(kFALSE),
59   fHeader(0),
60   fTracks(0),
61   fVertices(0),
62   fV0s(0),
63   fCascades(0),
64   fTracklets(0),
65   fJets(0),
66   fEmcalCells(0),
67   fPhosCells(0),
68   fCaloClusters(0),
69   fFmdClusters(0),
70   fPmdClusters(0),
71   fDimuons(0)
72 {
73   // default constructor
74 }
75
76 //______________________________________________________________________________
77 AliAODEvent::AliAODEvent(const AliAODEvent& aod):
78   AliVEvent(aod),
79   fAODObjects(new TList()),
80   fAODFolder(0),
81   fConnected(kFALSE),
82   fHeader(new AliAODHeader(*aod.fHeader)),
83   fTracks(new TClonesArray(*aod.fTracks)),
84   fVertices(new TClonesArray(*aod.fVertices)),
85   fV0s(new TClonesArray(*aod.fV0s)),
86   fCascades(new TClonesArray(*aod.fCascades)),
87   fTracklets(new AliAODTracklets(*aod.fTracklets)),
88   fJets(new TClonesArray(*aod.fJets)),
89   fEmcalCells(new AliAODCaloCells(*aod.fEmcalCells)),
90   fPhosCells(new AliAODCaloCells(*aod.fPhosCells)),
91   fCaloClusters(new TClonesArray(*aod.fCaloClusters)),
92   fFmdClusters(new TClonesArray(*aod.fFmdClusters)),
93   fPmdClusters(new TClonesArray(*aod.fPmdClusters)),
94   fDimuons(new TClonesArray(*aod.fDimuons))
95 {
96   // Copy constructor
97   AddObject(fHeader);
98   AddObject(fTracks);
99   AddObject(fVertices);
100   AddObject(fV0s);
101   AddObject(fCascades);
102   AddObject(fTracklets);
103   AddObject(fJets);
104   AddObject(fEmcalCells);
105   AddObject(fPhosCells);
106   AddObject(fCaloClusters);
107   AddObject(fFmdClusters);
108   AddObject(fPmdClusters);
109   AddObject(fDimuons);
110   fConnected = aod.fConnected;
111   GetStdContent();
112   CreateStdFolders();
113 }
114
115 //______________________________________________________________________________
116 AliAODEvent & AliAODEvent::operator=(const AliAODEvent& aod) {
117
118     // Assignment operator
119
120   if(&aod == this) return *this;
121   AliVEvent::operator=(aod);
122
123   // This assumes that the list is already created
124   // and that the virtual void Copy(Tobject&) function
125   // is correctly implemented in the derived class
126   // otherwise only TObject::Copy() will be used
127
128   if((fAODObjects->GetSize()==0)&&(aod.fAODObjects->GetSize()>=kAODListN)){
129     // We cover the case that we do not yet have the 
130     // standard content but the source has it
131     CreateStdContent();
132   }
133   
134   // Here we have the standard content without user additions, but the content is 
135   // not matching the aod source.
136   
137   // Iterate the list of source objects
138   TIter next(aod.GetList());
139   TObject *its = 0;
140   TString name;
141   while ((its = next())) {
142     name = its->GetName();
143     // Check if we have this object type in out list
144     TObject *mine = fAODObjects->FindObject(name);    
145     if(!mine) {
146       // We have to create the same type of object.
147       TClass* pClass=TClass::GetClass(its->ClassName());     
148       if (!pClass) {
149         AliWarning(Form("Can not find class description for entry %s (%s)\n",
150                    its->ClassName(), name.Data()));
151         continue;
152       }
153       mine=(TObject*)pClass->New();
154       if(!mine){
155         // not in this: can be added to list
156         AliWarning(Form("%s:%d Could not find %s for copying \n",
157                    (char*)__FILE__,__LINE__,name.Data()));
158         continue;
159       }  
160       if(mine->InheritsFrom("TNamed"))  {
161         ((TNamed*)mine)->SetName(name);
162       } else if(mine->InheritsFrom("TCollection")){
163         if(mine->InheritsFrom("TClonesArray")) {
164           TClonesArray *itscl = dynamic_cast<TClonesArray*>(its);
165           if (!itscl) {
166             AliWarning(Form("Class description for entry %s (%s) not TClonesArray\n",
167                    its->ClassName(), name.Data()));
168             continue;
169           
170           }
171                dynamic_cast<TClonesArray*>(mine)->SetClass(itscl->GetClass(), itscl->GetSize());
172         }
173         dynamic_cast<TCollection*>(mine)->SetName(name);
174       }
175       AliDebug(1, Form("adding object %s of type %s", mine->GetName(), mine->ClassName()));
176       AddObject(mine);
177     }
178     // Now we have an object of the same type and name, but different content.        
179     if(!its->InheritsFrom("TCollection")){
180       // simple objects (do they have a Copy method that calls operator= ?)
181       its->Copy(*mine);
182     } else if (its->InheritsFrom("TClonesArray")) {
183       // Create or expand the tclonesarray pointers
184       // so we can directly copy to the object
185       TClonesArray *its_tca = (TClonesArray*)its;
186       TClonesArray *mine_tca = (TClonesArray*)mine;
187       // this leaves the capacity of the TClonesArray the same
188       // except for a factor of 2 increase when size > capacity
189       // does not release any memory occupied by the tca
190       Int_t its_entries = its_tca->GetEntriesFast();
191       mine_tca->ExpandCreate(its_entries);
192       for(int i=0; i<its_entries; i++){
193         // copy 
194         TObject *mine_tca_obj = mine_tca->At(i);
195         TObject *its_tca_obj = its_tca->At(i);
196         // no need to delete first
197         // pointers within the class should be handled by Copy()...
198         // Can there be Empty slots?
199         its_tca_obj->Copy(*mine_tca_obj);
200       }
201     } else {
202       AliWarning(Form("%s:%d cannot copy TCollection \n",
203                       (char*)__FILE__,__LINE__));
204     }
205   }  
206   fConnected = aod.fConnected;
207   return *this;
208 }
209
210
211 //______________________________________________________________________________
212 AliAODEvent::~AliAODEvent() 
213 {
214 // destructor
215     if(fAODObjects&&!fConnected)
216     {
217         delete fAODObjects;
218         fAODObjects = 0;
219     }
220
221     delete fAODFolder;
222     fAODFolder = 0;
223 }
224
225 //______________________________________________________________________________
226 void AliAODEvent::AddObject(TObject* obj) 
227 {
228   // Add an object to the list of objects.
229   // Please be aware that in order to increase performance you should
230   // refrain from using TObjArrays (if possible). Use TClonesArrays, instead.
231   
232   fAODObjects->AddLast(obj);
233 }
234
235 //______________________________________________________________________________
236 void AliAODEvent::RemoveObject(TObject* obj) 
237 {
238   // Removes an object from the list of objects.
239   
240   fAODObjects->Remove(obj);
241 }
242
243 //______________________________________________________________________________
244 TObject *AliAODEvent::FindListObject(const char *objName) const
245 {
246   // Return the pointer to the object with the given name.
247
248   return fAODObjects->FindObject(objName);
249 }
250
251 //______________________________________________________________________________
252 void AliAODEvent::CreateStdContent() 
253 {
254   // create the standard AOD content and set pointers
255
256   // create standard objects and add them to the TList of objects
257   AddObject(new AliAODHeader());
258   AddObject(new TClonesArray("AliAODTrack", 0));
259   AddObject(new TClonesArray("AliAODVertex", 0));
260   AddObject(new TClonesArray("AliAODv0", 0));
261   AddObject(new TClonesArray("AliAODcascade", 0));
262   AddObject(new AliAODTracklets());
263   AddObject(new TClonesArray("AliAODJet", 0));
264   AddObject(new AliAODCaloCells());
265   AddObject(new AliAODCaloCells());
266   AddObject(new TClonesArray("AliAODCaloCluster", 0));
267   AddObject(new TClonesArray("AliAODFmdCluster", 0));
268   AddObject(new TClonesArray("AliAODPmdCluster", 0));
269   AddObject(new TClonesArray("AliAODDimuon", 0));
270   // set names
271   SetStdNames();
272
273   // read back pointers
274   GetStdContent();
275   CreateStdFolders();
276   return;
277 }
278
279 void  AliAODEvent::MakeEntriesReferencable()
280 {
281     // Make all entries referencable in a subsequent process
282     //
283     TIter next(fAODObjects);
284     TObject* obj;
285     while ((obj = next()))
286     {
287         if(obj->InheritsFrom("TCollection"))
288             {
289                 AssignIDtoCollection((TCollection*)obj);
290             }
291     }
292 }
293
294 //______________________________________________________________________________
295 void AliAODEvent::SetStdNames()
296 {
297   // introduce the standard naming
298
299   if(fAODObjects->GetEntries()==kAODListN){
300     for(int i = 0;i < fAODObjects->GetEntries();i++){
301       TObject *fObj = fAODObjects->At(i);
302       if(fObj->InheritsFrom("TNamed")){
303         ((TNamed*)fObj)->SetName(fAODListName[i]);
304       }
305       else if(fObj->InheritsFrom("TClonesArray")){
306         ((TClonesArray*)fObj)->SetName(fAODListName[i]);
307       }
308     }
309   }
310   else{
311     printf("%s:%d SetStdNames() Wrong number of Std Entries \n",(char*)__FILE__,__LINE__);
312   }
313
314
315 void AliAODEvent::CreateStdFolders()
316 {
317     // Create the standard folder structure
318   if(fAODFolder)delete fAODFolder;
319     fAODFolder = gROOT->GetRootFolder()->AddFolder("AOD", "AOD");
320     if(fAODObjects->GetEntries()==kAODListN){
321         for(int i = 0;i < fAODObjects->GetEntries();i++){
322             TObject *fObj = fAODObjects->At(i);
323             if(fObj->InheritsFrom("TClonesArray")){
324                 fAODFolder->AddFolder(fAODListName[i], fAODListName[i], (TCollection*) fObj);
325             } else {
326                 fAODFolder->AddFolder(fAODListName[i], fAODListName[i], 0);
327             }
328         }
329     }
330     else{
331         printf("%s:%d CreateStdFolders() Wrong number of Std Entries \n",(char*)__FILE__,__LINE__);
332     }
333
334
335 //______________________________________________________________________________
336 void AliAODEvent::GetStdContent()
337 {
338   // set pointers for standard content
339
340   fHeader        = (AliAODHeader*)fAODObjects->FindObject("header");
341   fTracks        = (TClonesArray*)fAODObjects->FindObject("tracks");
342   fVertices      = (TClonesArray*)fAODObjects->FindObject("vertices");
343   fV0s           = (TClonesArray*)fAODObjects->FindObject("v0s");
344   fCascades      = (TClonesArray*)fAODObjects->FindObject("cascades");
345   fTracklets     = (AliAODTracklets*)fAODObjects->FindObject("tracklets");
346   fJets          = (TClonesArray*)fAODObjects->FindObject("jets");
347   fEmcalCells    = (AliAODCaloCells*)fAODObjects->FindObject("emcalCells");
348   fPhosCells     = (AliAODCaloCells*)fAODObjects->FindObject("phosCells");
349   fCaloClusters  = (TClonesArray*)fAODObjects->FindObject("caloClusters");
350   fFmdClusters   = (TClonesArray*)fAODObjects->FindObject("fmdClusters");
351   fPmdClusters   = (TClonesArray*)fAODObjects->FindObject("pmdClusters");
352   fDimuons       = (TClonesArray*)fAODObjects->FindObject("dimuons");
353 }
354
355 //______________________________________________________________________________
356 void AliAODEvent::ResetStd(Int_t trkArrSize, 
357                            Int_t vtxArrSize, 
358                            Int_t v0ArrSize,
359                            Int_t cascadeArrSize,
360                            Int_t jetSize, 
361                            Int_t caloClusSize, 
362                            Int_t fmdClusSize, 
363                            Int_t pmdClusSize,
364                            Int_t dimuonArrSize
365                            )
366 {
367   // deletes content of standard arrays and resets size 
368   
369   fTracks->Delete();
370   if (trkArrSize > fTracks->GetSize()) 
371     fTracks->Expand(trkArrSize);
372
373   fVertices->Delete();
374   if (vtxArrSize > fVertices->GetSize()) 
375     fVertices->Expand(vtxArrSize);
376
377   fV0s->Delete();
378   if (v0ArrSize > fV0s->GetSize()) 
379     fV0s->Expand(v0ArrSize);
380   
381   fCascades->Delete();
382   if (cascadeArrSize > fCascades->GetSize()) 
383     fCascades->Expand(cascadeArrSize);
384   
385   fJets->Delete();
386   if (jetSize > fJets->GetSize())
387     fJets->Expand(jetSize);
388
389   fCaloClusters->Delete();
390   if (caloClusSize > fCaloClusters->GetSize()) 
391     fCaloClusters->Expand(caloClusSize);
392
393   fFmdClusters->Delete();
394   if (fmdClusSize > fFmdClusters->GetSize()) 
395     fFmdClusters->Expand(fmdClusSize);
396
397   fPmdClusters->Delete();
398   if (pmdClusSize > fPmdClusters->GetSize()) 
399     fPmdClusters->Expand(pmdClusSize);
400     
401   fDimuons->Delete();
402   if (dimuonArrSize > fDimuons->GetSize()) 
403     fDimuons->Expand(dimuonArrSize);
404
405   // Reset the tracklets
406   fTracklets->DeleteContainer();
407   fPhosCells->DeleteContainer();  
408   fEmcalCells->DeleteContainer();
409
410 }
411
412 void AliAODEvent::ClearStd()
413 {
414   // clears the standard arrays
415   fHeader        ->RemoveQTheta();
416   fTracks        ->Delete();
417   fVertices      ->Delete();
418   fV0s           ->Delete();
419   fCascades      ->Delete();
420   fTracklets     ->DeleteContainer();
421   fJets          ->Delete();
422   fEmcalCells    ->DeleteContainer();
423   fPhosCells     ->DeleteContainer();
424   fCaloClusters  ->Delete();
425   fFmdClusters   ->Clear();
426   fPmdClusters   ->Clear();
427   fDimuons       ->Clear();
428 }
429
430 //_________________________________________________________________
431 Int_t AliAODEvent::GetPHOSClusters(TRefArray *clusters) const
432 {
433   // fills the provided TRefArray with all found phos clusters
434   
435   clusters->Clear();
436   
437   AliAODCaloCluster *cl = 0;
438   Bool_t first = kTRUE;
439   for (Int_t i = 0; i < GetNCaloClusters() ; i++) {
440     if ( (cl = GetCaloCluster(i)) ) {
441       if (cl->IsPHOSCluster()){
442         if(first) {
443           new (clusters) TRefArray(TProcessID::GetProcessWithUID(cl)); 
444           first=kFALSE;
445         }
446         clusters->Add(cl);
447         //printf("IsPHOS cluster %d, E %2.3f Size: %d \n",i,cl->E(),clusters->GetEntriesFast());
448       }
449     }
450   }
451   return clusters->GetEntriesFast();
452 }
453
454 //_________________________________________________________________
455 Int_t AliAODEvent::GetEMCALClusters(TRefArray *clusters) const
456 {
457   // fills the provided TRefArray with all found emcal clusters
458
459   clusters->Clear();
460   AliAODCaloCluster *cl = 0;
461   Bool_t first = kTRUE;
462   for (Int_t i = 0; i < GetNCaloClusters(); i++) {
463     if ( (cl = GetCaloCluster(i)) ) {
464       if (cl->IsEMCALCluster()){
465         if(first) {
466           new (clusters) TRefArray(TProcessID::GetProcessWithUID(cl)); 
467           first=kFALSE;
468         }
469         clusters->Add(cl);
470         //printf("IsEMCal cluster %d, E %2.3f Size: %d \n",i,cl->E(),clusters->GetEntriesFast());
471       }
472     }
473   }
474   return clusters->GetEntriesFast();
475 }
476
477
478 //______________________________________________________________________________
479 Int_t AliAODEvent::GetMuonTracks(TRefArray *muonTracks) const
480 {
481   // fills the provided TRefArray with all found muon tracks
482
483   muonTracks->Clear();
484
485   AliAODTrack *track = 0;
486   for (Int_t iTrack = 0; iTrack < GetNTracks(); iTrack++) {
487     if ((track = GetTrack(iTrack))->IsMuonTrack()) {
488       muonTracks->Add(track);
489     }
490   }
491   
492   return muonTracks->GetEntriesFast();
493 }
494
495
496 //______________________________________________________________________________
497 Int_t AliAODEvent::GetNumberOfMuonTracks() const
498 {
499   // get number of muon tracks
500   Int_t nMuonTracks=0;
501   AliAODTrack *track = 0;
502   for (Int_t iTrack = 0; iTrack < GetNTracks(); iTrack++) {
503     if ((track = GetTrack(iTrack))->IsMuonTrack()) {
504        nMuonTracks++;
505     }
506   }
507   
508   return nMuonTracks;
509 }
510
511 void AliAODEvent::ReadFromTree(TTree *tree, Option_t* opt /*= ""*/)
512 {
513   // Connects aod event to tree
514   
515   if(!tree){
516     Printf("%s %d AliAODEvent::ReadFromTree() Zero Pointer to Tree \n",(char*)__FILE__,__LINE__);
517     return;
518   }
519   // load the TTree
520   if(!tree->GetTree())tree->LoadTree(0);
521
522   // Try to find AliAODEvent
523   AliAODEvent *aodEvent = 0;
524   aodEvent = (AliAODEvent*)tree->GetTree()->GetUserInfo()->FindObject("AliAODEvent");
525   if(aodEvent){
526     // Check if already connected to tree
527     TList* connectedList = (TList*) (tree->GetUserInfo()->FindObject("AODObjectsConnectedToTree"));
528     if (connectedList && (strcmp(opt, "reconnect"))) {
529         // If connected use the connected list of objects
530         fAODObjects->Delete();
531         fAODObjects = connectedList;
532         GetStdContent(); 
533         fConnected = kTRUE;
534         return;
535     } 
536     // Connect to tree
537     // prevent a memory leak when reading back the TList
538     // if (!(strcmp(opt, "reconnect"))) fAODObjects->Delete();
539
540     // create a new TList from the UserInfo TList... 
541     // copy constructor does not work...
542     fAODObjects = (TList*)(aodEvent->GetList()->Clone());
543     fAODObjects->SetOwner(kTRUE);
544     if(fAODObjects->GetEntries()<kAODListN){
545       printf("%s %d AliAODEvent::ReadFromTree() TList contains less than the standard contents %d < %d \n",
546              (char*)__FILE__,__LINE__,fAODObjects->GetEntries(),kAODListN);
547     }
548     //
549     // Let's find out whether we have friends
550     TList* friendL = tree->GetTree()->GetListOfFriends();
551     if (friendL) 
552     {
553         TIter next(friendL);
554         TFriendElement* fe;
555         while ((fe = (TFriendElement*)next())){
556             aodEvent = (AliAODEvent*)(fe->GetTree()->GetUserInfo()->FindObject("AliAODEvent"));
557             if (!aodEvent) {
558                 printf("No UserInfo on tree \n");
559             } else {
560
561                 TList* objL = (TList*)(aodEvent->GetList()->Clone());
562                 printf("Get list of object from tree %d !!\n", objL->GetEntries());
563                 TIter nextobject(objL);
564                 TObject* obj =  0;
565                 while((obj = nextobject()))
566                 {
567                     printf("Adding object from friend %s !\n", obj->GetName());
568                     fAODObjects->Add(obj);
569                 } // object "branch" loop
570             } // has userinfo  
571         } // friend loop
572     } // has friends    
573             
574
575 // set the branch addresses
576     TIter next(fAODObjects);
577     TNamed *el;
578     while((el=(TNamed*)next())){
579       TString bname(el->GetName());
580       // check if branch exists under this Name
581       TBranch *br = tree->GetTree()->GetBranch(bname.Data());
582       if(br){
583         tree->SetBranchAddress(bname.Data(),fAODObjects->GetObjectRef(el));
584       } else {
585           br = tree->GetBranch(Form("%s.",bname.Data()));
586           if(br){
587               tree->SetBranchAddress(Form("%s.",bname.Data()),fAODObjects->GetObjectRef(el));
588           }
589           else{
590               printf("%s %d AliAODEvent::ReadFromTree() No Branch found with Name %s. \n",
591                      (char*)__FILE__,__LINE__,bname.Data());
592           }     
593       }
594     }
595     GetStdContent();
596     // when reading back we are not owner of the list 
597     // must not delete it
598     fAODObjects->SetOwner(kTRUE);
599     fAODObjects->SetName("AODObjectsConnectedToTree");
600     // we are not owner of the list objects 
601     // must not delete it
602     tree->GetUserInfo()->Add(fAODObjects);
603     fConnected = kTRUE;
604   }// no aodEvent
605   else {
606     // we can't get the list from the user data, create standard content
607     // and set it by hand
608     CreateStdContent();
609     TIter next(fAODObjects);
610     TNamed *el;
611     while((el=(TNamed*)next())){
612       TString bname(el->GetName());    
613       tree->SetBranchAddress(bname.Data(),fAODObjects->GetObjectRef(el));
614     }
615     GetStdContent();
616     // when reading back we are not owner of the list 
617     // must not delete it
618     fAODObjects->SetOwner(kTRUE);
619   }
620 }
621 //______________________________________________________________________________
622 Int_t  AliAODEvent::GetNumberOfPileupVerticesSPD() const{
623   // count number of SPD pileup vertices
624   Int_t nVertices=GetNumberOfVertices();
625   Int_t nPileupVertices=0;
626   for(Int_t iVert=0; iVert<nVertices; iVert++){
627     AliAODVertex *v=GetVertex(iVert);
628     if(v->GetType()==AliAODVertex::kPileupSPD) nPileupVertices++;
629   }
630   return nPileupVertices;
631 }
632 //______________________________________________________________________________
633 Int_t  AliAODEvent::GetNumberOfPileupVerticesTracks() const{
634   // count number of track pileup vertices
635   Int_t nVertices=GetNumberOfVertices();
636   Int_t nPileupVertices=0;
637   for(Int_t iVert=0; iVert<nVertices; iVert++){
638     AliAODVertex *v=GetVertex(iVert);
639     if(v->GetType()==AliAODVertex::kPileupTracks) nPileupVertices++;
640   }
641   return nPileupVertices;
642 }
643 //______________________________________________________________________________
644 AliAODVertex* AliAODEvent::GetPrimaryVertexSPD() const{
645   //
646   Int_t nVertices=GetNumberOfVertices();
647   for(Int_t iVert=0; iVert<nVertices; iVert++){
648     AliAODVertex *v=GetVertex(iVert);
649     if(v->GetType()==AliAODVertex::kMainSPD) return v;
650   }
651   return 0;
652 }
653 //______________________________________________________________________________
654 AliAODVertex* AliAODEvent::GetPileupVertexSPD(Int_t iV) const{
655   //
656   Int_t nVertices=GetNumberOfVertices();
657   Int_t counter=0;
658   for(Int_t iVert=0; iVert<nVertices; iVert++){
659     AliAODVertex *v=GetVertex(iVert);
660     if(v->GetType()==AliAODVertex::kPileupSPD){
661       if(counter==iV) return v;
662       ++counter;
663     }
664   }
665   return 0;
666 }
667 //______________________________________________________________________________
668 AliAODVertex* AliAODEvent::GetPileupVertexTracks(Int_t iV) const{
669   //
670   Int_t nVertices=GetNumberOfVertices();
671   Int_t counter=0;
672   for(Int_t iVert=0; iVert<nVertices; iVert++){
673     AliAODVertex *v=GetVertex(iVert);
674     if(v->GetType()==AliAODVertex::kPileupTracks){
675       if(counter==iV) return v;
676       ++counter;
677     }
678   }
679   return 0;
680 }
681 //______________________________________________________________________________
682 Bool_t  AliAODEvent::IsPileupFromSPD(Int_t minContributors, 
683                                      Double_t minZdist, 
684                                      Double_t nSigmaZdist, 
685                                      Double_t nSigmaDiamXY, 
686                                      Double_t nSigmaDiamZ) const{
687   //
688   // This function checks if there was a pile up
689   // reconstructed with SPD
690   //
691   AliAODVertex *mainV=GetPrimaryVertexSPD();
692   if(!mainV) return kFALSE;
693   Int_t nc1=mainV->GetNContributors();
694   if(nc1<1) return kFALSE;
695   Int_t nPileVert=GetNumberOfPileupVerticesSPD();
696   if(nPileVert==0) return kFALSE;
697   Int_t nVertices=GetNumberOfVertices();
698   
699   for(Int_t iVert=0; iVert<nVertices; iVert++){
700     AliAODVertex *pv=GetVertex(iVert);
701     if(pv->GetType()!=AliAODVertex::kPileupSPD) continue;
702     Int_t nc2=pv->GetNContributors();
703     if(nc2>=minContributors){
704       Double_t z1=mainV->GetZ();
705       Double_t z2=pv->GetZ();
706       Double_t distZ=TMath::Abs(z2-z1);
707       Double_t distZdiam=TMath::Abs(z2-GetDiamondZ());
708       Double_t cutZdiam=nSigmaDiamZ*GetSigma2DiamondZ();
709       if(GetSigma2DiamondZ()<0.0001)cutZdiam=99999.; //protection for missing z diamond information
710       if(distZ>minZdist && distZdiam<cutZdiam){
711         Double_t x2=pv->GetX();
712         Double_t y2=pv->GetY();
713         Double_t distXdiam=TMath::Abs(x2-GetDiamondX());
714         Double_t distYdiam=TMath::Abs(y2-GetDiamondY());
715         Double_t cov1[6],cov2[6];       
716         mainV->GetCovarianceMatrix(cov1);
717         pv->GetCovarianceMatrix(cov2);
718         Double_t errxDist=TMath::Sqrt(cov2[0]+GetSigma2DiamondX());
719         Double_t erryDist=TMath::Sqrt(cov2[2]+GetSigma2DiamondY());
720         Double_t errzDist=TMath::Sqrt(cov1[5]+cov2[5]);
721         Double_t cutXdiam=nSigmaDiamXY*errxDist;
722         if(GetSigma2DiamondX()<0.0001)cutXdiam=99999.; //protection for missing diamond information
723         Double_t cutYdiam=nSigmaDiamXY*erryDist;
724         if(GetSigma2DiamondY()<0.0001)cutYdiam=99999.; //protection for missing diamond information
725         if( (distXdiam<cutXdiam) && (distYdiam<cutYdiam) && (distZ>nSigmaZdist*errzDist) ){
726           return kTRUE;
727         }
728       }
729     }
730   }
731   return kFALSE;
732 }
733
734 //______________________________________________________________________________
735 void AliAODEvent::Print(Option_t *) const
736 {
737   // Print the names of the all branches
738   TIter next(fAODObjects);
739   TNamed *el;
740   Printf(">>>>>  AOD  Content <<<<<");    
741   while((el=(TNamed*)next())){
742     Printf(">> %s ",el->GetName());      
743   }
744   Printf(">>>>>                <<<<<");    
745   
746   return;
747 }
748
749 void AliAODEvent::AssignIDtoCollection(TCollection* col)
750 {
751     // Static method which assigns a ID to each object in a collection
752     // In this way the objects are marked as referenced and written with 
753     // an ID. This has the advantage that TRefs to this objects can be 
754     // written by a subsequent process.
755     TIter next(col);
756     TObject* obj;
757     while ((obj = next()))
758         TProcessID::AssignID(obj);
759 }