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