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