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