]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliAODEvent.cxx
Fix for Coverity warning 16208
[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(fAODObjects&&!fConnected)
220     {
221         delete fAODObjects;
222         fAODObjects = 0;
223     }
224
225     delete fAODFolder;
226     fAODFolder = 0;
227 }
228
229 //______________________________________________________________________________
230 void AliAODEvent::AddObject(TObject* obj) 
231 {
232   // Add an object to the list of objects.
233   // Please be aware that in order to increase performance you should
234   // refrain from using TObjArrays (if possible). Use TClonesArrays, instead.
235   
236   if ( !fAODObjects->FindObject(obj) ) 
237   {
238     fAODObjects->AddLast(obj);
239   }
240 }
241
242 //______________________________________________________________________________
243 void AliAODEvent::RemoveObject(TObject* obj) 
244 {
245   // Removes an object from the list of objects.
246   
247   fAODObjects->Remove(obj);
248 }
249
250 //______________________________________________________________________________
251 TObject *AliAODEvent::FindListObject(const char *objName) const
252 {
253   // Return the pointer to the object with the given name.
254
255   return fAODObjects->FindObject(objName);
256 }
257
258 //______________________________________________________________________________
259 void AliAODEvent::CreateStdContent() 
260 {
261   // create the standard AOD content and set pointers
262
263   // create standard objects and add them to the TList of objects
264   AddObject(new AliAODHeader());
265   AddObject(new TClonesArray("AliAODTrack", 0));
266   AddObject(new TClonesArray("AliAODVertex", 0));
267   AddObject(new TClonesArray("AliAODv0", 0));
268   AddObject(new TClonesArray("AliAODcascade", 0));
269   AddObject(new AliAODTracklets());
270   AddObject(new TClonesArray("AliAODJet", 0));
271   AddObject(new AliAODCaloCells());
272   AddObject(new AliAODCaloCells());
273   AddObject(new TClonesArray("AliAODCaloCluster", 0));
274   AddObject(new TClonesArray("AliAODFmdCluster", 0));
275   AddObject(new TClonesArray("AliAODPmdCluster", 0));
276   AddObject(new TClonesArray("AliAODDimuon", 0));
277   AddObject(new AliAODVZERO());
278   // set names
279   SetStdNames();
280
281   // read back pointers
282   GetStdContent();
283   CreateStdFolders();
284   return;
285 }
286
287 void  AliAODEvent::MakeEntriesReferencable()
288 {
289     // Make all entries referencable in a subsequent process
290     //
291     TIter next(fAODObjects);
292     TObject* obj;
293     while ((obj = next()))
294     {
295         if(obj->InheritsFrom("TCollection"))
296             {
297                 AssignIDtoCollection((TCollection*)obj);
298             }
299     }
300 }
301
302 //______________________________________________________________________________
303 void AliAODEvent::SetStdNames()
304 {
305   // introduce the standard naming
306
307   if(fAODObjects->GetEntries()==kAODListN){
308     for(int i = 0;i < fAODObjects->GetEntries();i++){
309       TObject *fObj = fAODObjects->At(i);
310       if(fObj->InheritsFrom("TNamed")){
311         ((TNamed*)fObj)->SetName(fAODListName[i]);
312       }
313       else if(fObj->InheritsFrom("TClonesArray")){
314         ((TClonesArray*)fObj)->SetName(fAODListName[i]);
315       }
316     }
317   }
318   else{
319     printf("%s:%d SetStdNames() Wrong number of Std Entries \n",(char*)__FILE__,__LINE__);
320   }
321
322
323 void AliAODEvent::CreateStdFolders()
324 {
325     // Create the standard folder structure
326   if(fAODFolder)delete fAODFolder;
327     fAODFolder = gROOT->GetRootFolder()->AddFolder("AOD", "AOD");
328     if(fAODObjects->GetEntries()==kAODListN){
329         for(int i = 0;i < fAODObjects->GetEntries();i++){
330             TObject *fObj = fAODObjects->At(i);
331             if(fObj->InheritsFrom("TClonesArray")){
332                 fAODFolder->AddFolder(fAODListName[i], fAODListName[i], (TCollection*) fObj);
333             } else {
334                 fAODFolder->AddFolder(fAODListName[i], fAODListName[i], 0);
335             }
336         }
337     }
338     else{
339         printf("%s:%d CreateStdFolders() Wrong number of Std Entries \n",(char*)__FILE__,__LINE__);
340     }
341
342
343 //______________________________________________________________________________
344 void AliAODEvent::GetStdContent()
345 {
346   // set pointers for standard content
347
348   fHeader        = (AliAODHeader*)fAODObjects->FindObject("header");
349   fTracks        = (TClonesArray*)fAODObjects->FindObject("tracks");
350   fVertices      = (TClonesArray*)fAODObjects->FindObject("vertices");
351   fV0s           = (TClonesArray*)fAODObjects->FindObject("v0s");
352   fCascades      = (TClonesArray*)fAODObjects->FindObject("cascades");
353   fTracklets     = (AliAODTracklets*)fAODObjects->FindObject("tracklets");
354   fJets          = (TClonesArray*)fAODObjects->FindObject("jets");
355   fEmcalCells    = (AliAODCaloCells*)fAODObjects->FindObject("emcalCells");
356   fPhosCells     = (AliAODCaloCells*)fAODObjects->FindObject("phosCells");
357   fCaloClusters  = (TClonesArray*)fAODObjects->FindObject("caloClusters");
358   fFmdClusters   = (TClonesArray*)fAODObjects->FindObject("fmdClusters");
359   fPmdClusters   = (TClonesArray*)fAODObjects->FindObject("pmdClusters");
360   fDimuons       = (TClonesArray*)fAODObjects->FindObject("dimuons");
361   fAODVZERO      = (AliAODVZERO*)fAODObjects->FindObject("AliAODVZERO");
362 }
363
364 //______________________________________________________________________________
365 void AliAODEvent::ResetStd(Int_t trkArrSize, 
366                            Int_t vtxArrSize, 
367                            Int_t v0ArrSize,
368                            Int_t cascadeArrSize,
369                            Int_t jetSize, 
370                            Int_t caloClusSize, 
371                            Int_t fmdClusSize, 
372                            Int_t pmdClusSize,
373                            Int_t dimuonArrSize
374                            )
375 {
376   // deletes content of standard arrays and resets size 
377   
378   if (fTracks) {
379     fTracks->Delete();
380     if (trkArrSize > fTracks->GetSize()) 
381       fTracks->Expand(trkArrSize);
382   }
383   if (fVertices) {
384     fVertices->Delete();
385     if (vtxArrSize > fVertices->GetSize()) 
386       fVertices->Expand(vtxArrSize);
387   }
388   if (fV0s) {
389     fV0s->Delete();
390     if (v0ArrSize > fV0s->GetSize()) 
391       fV0s->Expand(v0ArrSize);
392   }
393   if (fCascades) {
394     fCascades->Delete();
395     if (cascadeArrSize > fCascades->GetSize()) 
396       fCascades->Expand(cascadeArrSize);
397   }
398   if (fJets) {
399     fJets->Delete();
400     if (jetSize > fJets->GetSize())
401       fJets->Expand(jetSize);
402   }
403   if (fCaloClusters) {
404     fCaloClusters->Delete();
405     if (caloClusSize > fCaloClusters->GetSize()) 
406       fCaloClusters->Expand(caloClusSize);
407   }
408   if (fFmdClusters) {
409     fFmdClusters->Delete();
410     if (fmdClusSize > fFmdClusters->GetSize()) 
411       fFmdClusters->Expand(fmdClusSize);
412   }
413   if (fPmdClusters) {
414     fPmdClusters->Delete();
415     if (pmdClusSize > fPmdClusters->GetSize()) 
416       fPmdClusters->Expand(pmdClusSize);
417   }
418   if (fDimuons) {
419     fDimuons->Delete();
420     if (dimuonArrSize > fDimuons->GetSize()) 
421       fDimuons->Expand(dimuonArrSize);
422   }
423   if (fTracklets)
424     fTracklets->DeleteContainer();
425   if (fPhosCells)
426     fPhosCells->DeleteContainer();  
427   if (fEmcalCells)
428     fEmcalCells->DeleteContainer();
429 }
430
431 void AliAODEvent::ClearStd()
432 {
433   // clears the standard arrays
434   if (fHeader)
435     fHeader        ->Clear();
436   if (fTracks)
437     fTracks        ->Delete();
438   if (fVertices)
439     fVertices      ->Delete();
440   if (fV0s)
441     fV0s           ->Delete();
442   if (fCascades)
443     fCascades      ->Delete();
444   if (fTracklets)
445     fTracklets     ->DeleteContainer();
446   if (fJets)
447     fJets          ->Delete();
448   if (fEmcalCells)
449     fEmcalCells    ->DeleteContainer();
450   if (fPhosCells)
451     fPhosCells     ->DeleteContainer();
452   if (fCaloClusters)
453     fCaloClusters  ->Delete();
454   if (fFmdClusters)
455     fFmdClusters   ->Clear();
456   if (fPmdClusters)
457     fPmdClusters   ->Clear();
458   if (fDimuons)
459     fDimuons       ->Clear();
460 }
461
462 //_________________________________________________________________
463 Int_t AliAODEvent::GetPHOSClusters(TRefArray *clusters) const
464 {
465   // fills the provided TRefArray with all found phos clusters
466   
467   clusters->Clear();
468   
469   AliAODCaloCluster *cl = 0;
470   Bool_t first = kTRUE;
471   for (Int_t i = 0; i < GetNumberOfCaloClusters() ; i++) {
472     if ( (cl = GetCaloCluster(i)) ) {
473       if (cl->IsPHOS()){
474         if(first) {
475           new (clusters) TRefArray(TProcessID::GetProcessWithUID(cl)); 
476           first=kFALSE;
477         }
478         clusters->Add(cl);
479         //printf("IsPHOS cluster %d, E %2.3f Size: %d \n",i,cl->E(),clusters->GetEntriesFast());
480       }
481     }
482   }
483   return clusters->GetEntriesFast();
484 }
485
486 //_________________________________________________________________
487 Int_t AliAODEvent::GetEMCALClusters(TRefArray *clusters) const
488 {
489   // fills the provided TRefArray with all found emcal clusters
490
491   clusters->Clear();
492   AliAODCaloCluster *cl = 0;
493   Bool_t first = kTRUE;
494   for (Int_t i = 0; i < GetNumberOfCaloClusters(); i++) {
495     if ( (cl = GetCaloCluster(i)) ) {
496       if (cl->IsEMCAL()){
497         if(first) {
498           new (clusters) TRefArray(TProcessID::GetProcessWithUID(cl)); 
499           first=kFALSE;
500         }
501         clusters->Add(cl);
502         //printf("IsEMCal cluster %d, E %2.3f Size: %d \n",i,cl->E(),clusters->GetEntriesFast());
503       }
504     }
505   }
506   return clusters->GetEntriesFast();
507 }
508
509
510 //______________________________________________________________________________
511 Int_t AliAODEvent::GetMuonTracks(TRefArray *muonTracks) const
512 {
513   // fills the provided TRefArray with all found muon tracks
514
515   muonTracks->Clear();
516
517   AliAODTrack *track = 0;
518   for (Int_t iTrack = 0; iTrack < GetNTracks(); iTrack++) {
519     track = GetTrack(iTrack);
520     if (track->IsMuonTrack()) {
521       muonTracks->Add(track);
522     }
523   }
524   
525   return muonTracks->GetEntriesFast();
526 }
527
528
529 //______________________________________________________________________________
530 Int_t AliAODEvent::GetNumberOfMuonTracks() const
531 {
532   // get number of muon tracks
533   Int_t nMuonTracks=0;
534   for (Int_t iTrack = 0; iTrack < GetNTracks(); iTrack++) {
535     if ((GetTrack(iTrack))->IsMuonTrack()) {
536        nMuonTracks++;
537     }
538   }
539   
540   return nMuonTracks;
541 }
542
543 void AliAODEvent::ReadFromTree(TTree *tree, Option_t* opt /*= ""*/)
544 {
545     // Connects aod event to tree
546   
547   if(!tree){
548     Printf("%s %d AliAODEvent::ReadFromTree() Zero Pointer to Tree \n",(char*)__FILE__,__LINE__);
549     return;
550   }
551     // load the TTree
552   if(!tree->GetTree())tree->LoadTree(0);
553   
554     // Try to find AliAODEvent
555   AliAODEvent *aodEvent = 0;
556   aodEvent = (AliAODEvent*)tree->GetTree()->GetUserInfo()->FindObject("AliAODEvent");
557   if(aodEvent){
558       // Check if already connected to tree
559     TList* connectedList = (TList*) (tree->GetUserInfo()->FindObject("AODObjectsConnectedToTree"));
560     if (connectedList && (strcmp(opt, "reconnect"))) {
561         // If connected use the connected list of objects
562         fAODObjects->Delete();
563         fAODObjects = connectedList;
564         GetStdContent(); 
565         fConnected = kTRUE;
566         return;
567     } 
568       // Connect to tree
569       // prevent a memory leak when reading back the TList
570 //      if (!(strcmp(opt, "reconnect"))) fAODObjects->Delete();
571     
572       // create a new TList from the UserInfo TList... 
573       // copy constructor does not work...
574     fAODObjects = (TList*)(aodEvent->GetList()->Clone());
575     fAODObjects->SetOwner(kTRUE);
576     if(fAODObjects->GetEntries()<kAODListN)
577     {
578       AliWarning(Form("AliAODEvent::ReadFromTree() TList contains less than the standard contents %d < %d"
579                       " That might be fine though (at least for filtered AODs)",fAODObjects->GetEntries(),kAODListN));
580     }
581       //
582       // Let's find out whether we have friends
583     TList* friendL = tree->GetTree()->GetListOfFriends();
584     if (friendL) 
585     {
586       TIter next(friendL);
587       TFriendElement* fe;
588       while ((fe = (TFriendElement*)next())){
589         aodEvent = (AliAODEvent*)(fe->GetTree()->GetUserInfo()->FindObject("AliAODEvent"));
590         if (!aodEvent) {
591           printf("No UserInfo on tree \n");
592         } else {
593           
594           TList* objL = (TList*)(aodEvent->GetList()->Clone());
595           printf("Get list of object from tree %d !!\n", objL->GetEntries());
596           TIter nextobject(objL);
597           TObject* obj =  0;
598           while((obj = nextobject()))
599           {
600             printf("Adding object from friend %s !\n", obj->GetName());
601             fAODObjects->Add(obj);
602           } // object "branch" loop
603         } // has userinfo  
604       } // friend loop
605     } // has friends    
606       // set the branch addresses
607     TIter next(fAODObjects);
608     TNamed *el;
609     while((el=(TNamed*)next())){
610       TString bname(el->GetName());
611         // check if branch exists under this Name
612       TBranch *br = tree->GetTree()->GetBranch(bname.Data());
613       if(br){
614         tree->SetBranchAddress(bname.Data(),fAODObjects->GetObjectRef(el));
615       } else {
616         br = tree->GetBranch(Form("%s.",bname.Data()));
617         if(br){
618           tree->SetBranchAddress(Form("%s.",bname.Data()),fAODObjects->GetObjectRef(el));
619         }
620         else{
621           printf("%s %d AliAODEvent::ReadFromTree() No Branch found with Name %s. \n",
622                  (char*)__FILE__,__LINE__,bname.Data());
623         }       
624       }
625     }
626     GetStdContent();
627       // when reading back we are not owner of the list 
628       // must not delete it
629     fAODObjects->SetOwner(kTRUE);
630     fAODObjects->SetName("AODObjectsConnectedToTree");
631       // we are not owner of the list objects 
632       // must not delete it
633     tree->GetUserInfo()->Add(fAODObjects);
634     fConnected = kTRUE;
635   }// no aodEvent
636   else {
637       // we can't get the list from the user data, create standard content
638       // and set it by hand
639     CreateStdContent();
640     TIter next(fAODObjects);
641     TNamed *el;
642     while((el=(TNamed*)next())){
643       TString bname(el->GetName());    
644       tree->SetBranchAddress(bname.Data(),fAODObjects->GetObjectRef(el));
645     }
646     GetStdContent();
647       // when reading back we are not owner of the list 
648       // must not delete it
649     fAODObjects->SetOwner(kTRUE);
650   }
651 }
652   //______________________________________________________________________________
653 Int_t  AliAODEvent::GetNumberOfPileupVerticesSPD() const{
654   // count number of SPD pileup vertices
655   Int_t nVertices=GetNumberOfVertices();
656   Int_t nPileupVertices=0;
657   for(Int_t iVert=0; iVert<nVertices; iVert++){
658     AliAODVertex *v=GetVertex(iVert);
659     if(v->GetType()==AliAODVertex::kPileupSPD) nPileupVertices++;
660   }
661   return nPileupVertices;
662 }
663 //______________________________________________________________________________
664 Int_t  AliAODEvent::GetNumberOfPileupVerticesTracks() const{
665   // count number of track pileup vertices
666   Int_t nVertices=GetNumberOfVertices();
667   Int_t nPileupVertices=0;
668   for(Int_t iVert=0; iVert<nVertices; iVert++){
669     AliAODVertex *v=GetVertex(iVert);
670     if(v->GetType()==AliAODVertex::kPileupTracks) nPileupVertices++;
671   }
672   return nPileupVertices;
673 }
674 //______________________________________________________________________________
675 AliAODVertex* AliAODEvent::GetPrimaryVertexSPD() const{
676   //
677   Int_t nVertices=GetNumberOfVertices();
678   for(Int_t iVert=0; iVert<nVertices; iVert++){
679     AliAODVertex *v=GetVertex(iVert);
680     if(v->GetType()==AliAODVertex::kMainSPD) return v;
681   }
682   return 0;
683 }
684 //______________________________________________________________________________
685 AliAODVertex* AliAODEvent::GetPileupVertexSPD(Int_t iV) const{
686   //
687   Int_t nVertices=GetNumberOfVertices();
688   Int_t counter=0;
689   for(Int_t iVert=0; iVert<nVertices; iVert++){
690     AliAODVertex *v=GetVertex(iVert);
691     if(v->GetType()==AliAODVertex::kPileupSPD){
692       if(counter==iV) return v;
693       ++counter;
694     }
695   }
696   return 0;
697 }
698 //______________________________________________________________________________
699 AliAODVertex* AliAODEvent::GetPileupVertexTracks(Int_t iV) const{
700   //
701   Int_t nVertices=GetNumberOfVertices();
702   Int_t counter=0;
703   for(Int_t iVert=0; iVert<nVertices; iVert++){
704     AliAODVertex *v=GetVertex(iVert);
705     if(v->GetType()==AliAODVertex::kPileupTracks){
706       if(counter==iV) return v;
707       ++counter;
708     }
709   }
710   return 0;
711 }
712 //______________________________________________________________________________
713 Bool_t  AliAODEvent::IsPileupFromSPD(Int_t minContributors, 
714                                      Double_t minZdist, 
715                                      Double_t nSigmaZdist, 
716                                      Double_t nSigmaDiamXY, 
717                                      Double_t nSigmaDiamZ) const{
718   //
719   // This function checks if there was a pile up
720   // reconstructed with SPD
721   //
722   AliAODVertex *mainV=GetPrimaryVertexSPD();
723   if(!mainV) return kFALSE;
724   Int_t nc1=mainV->GetNContributors();
725   if(nc1<1) return kFALSE;
726   Int_t nPileVert=GetNumberOfPileupVerticesSPD();
727   if(nPileVert==0) return kFALSE;
728   Int_t nVertices=GetNumberOfVertices();
729   
730   for(Int_t iVert=0; iVert<nVertices; iVert++){
731     AliAODVertex *pv=GetVertex(iVert);
732     if(pv->GetType()!=AliAODVertex::kPileupSPD) continue;
733     Int_t nc2=pv->GetNContributors();
734     if(nc2>=minContributors){
735       Double_t z1=mainV->GetZ();
736       Double_t z2=pv->GetZ();
737       Double_t distZ=TMath::Abs(z2-z1);
738       Double_t distZdiam=TMath::Abs(z2-GetDiamondZ());
739       Double_t cutZdiam=nSigmaDiamZ*TMath::Sqrt(GetSigma2DiamondZ());
740       if(GetSigma2DiamondZ()<0.0001)cutZdiam=99999.; //protection for missing z diamond information
741       if(distZ>minZdist && distZdiam<cutZdiam){
742         Double_t x2=pv->GetX();
743         Double_t y2=pv->GetY();
744         Double_t distXdiam=TMath::Abs(x2-GetDiamondX());
745         Double_t distYdiam=TMath::Abs(y2-GetDiamondY());
746         Double_t cov1[6],cov2[6];       
747         mainV->GetCovarianceMatrix(cov1);
748         pv->GetCovarianceMatrix(cov2);
749         Double_t errxDist=TMath::Sqrt(cov2[0]+GetSigma2DiamondX());
750         Double_t erryDist=TMath::Sqrt(cov2[2]+GetSigma2DiamondY());
751         Double_t errzDist=TMath::Sqrt(cov1[5]+cov2[5]);
752         Double_t cutXdiam=nSigmaDiamXY*errxDist;
753         if(GetSigma2DiamondX()<0.0001)cutXdiam=99999.; //protection for missing diamond information
754         Double_t cutYdiam=nSigmaDiamXY*erryDist;
755         if(GetSigma2DiamondY()<0.0001)cutYdiam=99999.; //protection for missing diamond information
756         if( (distXdiam<cutXdiam) && (distYdiam<cutYdiam) && (distZ>nSigmaZdist*errzDist) ){
757           return kTRUE;
758         }
759       }
760     }
761   }
762   return kFALSE;
763 }
764
765 //______________________________________________________________________________
766 void AliAODEvent::Print(Option_t *) const
767 {
768   // Print the names of the all branches
769   TIter next(fAODObjects);
770   TNamed *el;
771   Printf(">>>>>  AOD  Content <<<<<");    
772   while((el=(TNamed*)next())){
773     Printf(">> %s ",el->GetName());      
774   }
775   Printf(">>>>>                <<<<<");    
776   
777   return;
778 }
779
780 void AliAODEvent::AssignIDtoCollection(TCollection* col)
781 {
782     // Static method which assigns a ID to each object in a collection
783     // In this way the objects are marked as referenced and written with 
784     // an ID. This has the advantage that TRefs to this objects can be 
785     // written by a subsequent process.
786     TIter next(col);
787     TObject* obj;
788     while ((obj = next()))
789         TProcessID::AssignID(obj);
790 }
791
792 Bool_t AliAODEvent::IsPileupFromSPDInMultBins() const {
793     Int_t nTracklets=GetTracklets()->GetNumberOfTracklets();
794     if(nTracklets<20) return IsPileupFromSPD(3,0.8);
795     else if(nTracklets<50) return IsPileupFromSPD(4,0.8);
796     else return IsPileupFromSPD(5,0.8);
797 }
798