]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliESDEvent.cxx
New macro to keep track of timing performances of the segmentation methods (Laurent)
[u/mrichter/AliRoot.git] / STEER / AliESDEvent.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, 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 //           Implementation of the AliESDEvent class
20 //   This is the class to deal with during the physics analysis of data.
21 //   It also ensures the backward compatibility with the old ESD format.
22 /*
23    AliESDEvent *ev= new AliESDEvent();
24    ev->ReadFromTree(esdTree);
25    ...
26     for (Int_t i=0; i<nev; i++) {
27       esdTree->GetEntry(i);
28       if(ev->GetAliESDOld())ev->CopyFromOldESD();
29 */
30 //   The AliESDInputHandler does this automatically for you
31 //
32 // Origin: Christian Klein-Boesing, CERN, Christian.Klein-Boesing@cern.ch
33 //-----------------------------------------------------------------
34
35 #include "TList.h"
36 #include "TRefArray.h"
37 #include <TNamed.h>
38
39 #include "AliESDEvent.h"
40 #include "AliESDfriend.h"
41 #include "AliESDVZERO.h"
42 #include "AliESDFMD.h"
43 #include "AliESD.h"
44 #include "AliESDMuonTrack.h"
45 #include "AliESDPmdTrack.h"
46 #include "AliESDTrdTrack.h"
47 #include "AliESDVertex.h"
48 #include "AliESDcascade.h"
49 #include "AliESDPmdTrack.h"
50 #include "AliESDTrdTrack.h"
51 #include "AliESDVertex.h"
52 #include "AliESDcascade.h"
53 #include "AliESDkink.h"
54 #include "AliESDtrack.h"
55 #include "AliESDHLTtrack.h"
56 #include "AliESDCaloCluster.h"
57 #include "AliESDCaloCells.h"
58 #include "AliESDv0.h"
59 #include "AliESDFMD.h"
60 #include "AliESDVZERO.h"
61 #include "AliMultiplicity.h"
62 #include "AliRawDataErrorLog.h"
63 #include "AliLog.h"
64 #include "AliESDACORDE.h"
65 ClassImp(AliESDEvent)
66
67
68
69 // here we define the names, some classes are no TNamed, therefore the classnames 
70 // are the Names
71   const char* AliESDEvent::fgkESDListName[kESDListN] = {"AliESDRun",
72                                                        "AliESDHeader",
73                                                        "AliESDZDC",
74                                                        "AliESDFMD",
75                                                        "AliESDVZERO",
76                                                        "AliESDTZERO",
77                                                        "TPCVertex",
78                                                        "SPDVertex",
79                                                        "PrimaryVertex",
80                                                        "AliMultiplicity",
81                                                        "PHOSTrigger",
82                                                        "EMCALTrigger",
83                                                        "Tracks",
84                                                        "MuonTracks",
85                                                        "PmdTracks",
86                                                        "TrdTracks",
87                                                        "V0s",
88                                                        "Cascades",
89                                                        "Kinks",
90                                                        "CaloClusters",
91                                                       "EMCALCells",
92                                                       "PHOSCells",
93                                                        "AliRawDataErrorLogs",
94                                                        "AliESDACORDE"};
95
96 //______________________________________________________________________________
97 AliESDEvent::AliESDEvent():
98   AliVEvent(),
99   fESDObjects(new TList()),
100   fESDRun(0),
101   fHeader(0),
102   fESDZDC(0),
103   fESDFMD(0),
104   fESDVZERO(0),
105   fESDTZERO(0),
106   fTPCVertex(0),
107   fSPDVertex(0),
108   fPrimaryVertex(0),
109   fSPDMult(0),
110   fPHOSTrigger(0),
111   fEMCALTrigger(0),
112   fESDACORDE(0),
113   fTracks(0),
114   fMuonTracks(0),
115   fPmdTracks(0),
116   fTrdTracks(0),
117   fV0s(0),  
118   fCascades(0),
119   fKinks(0),
120   fCaloClusters(0),
121   fEMCALCells(0), fPHOSCells(0),
122   fErrorLogs(0),
123   fESDOld(0),
124   fESDFriendOld(0),
125   fConnected(kFALSE),
126   fUseOwnList(kFALSE),
127   fEMCALClusters(0), 
128   fFirstEMCALCluster(-1),
129   fPHOSClusters(0), 
130   fFirstPHOSCluster(-1)
131 {
132 }
133 //______________________________________________________________________________
134 AliESDEvent::AliESDEvent(const AliESDEvent& esd):
135   AliVEvent(esd),
136   fESDObjects(new TList()),
137   fESDRun(new AliESDRun(*esd.fESDRun)),
138   fHeader(new AliESDHeader(*esd.fHeader)),
139   fESDZDC(new AliESDZDC(*esd.fESDZDC)),
140   fESDFMD(new AliESDFMD(*esd.fESDFMD)),
141   fESDVZERO(new AliESDVZERO(*esd.fESDVZERO)),
142   fESDTZERO(new AliESDTZERO(*esd.fESDTZERO)),
143   fTPCVertex(new AliESDVertex(*esd.fTPCVertex)),
144   fSPDVertex(new AliESDVertex(*esd.fSPDVertex)),
145   fPrimaryVertex(new AliESDVertex(*esd.fPrimaryVertex)),
146   fSPDMult(new AliMultiplicity(*esd.fSPDMult)),
147   fPHOSTrigger(new AliESDCaloTrigger(*esd.fPHOSTrigger)),
148   fEMCALTrigger(new AliESDCaloTrigger(*esd.fEMCALTrigger)),
149   fESDACORDE(new AliESDACORDE(*esd.fESDACORDE)),
150   fTracks(new TClonesArray(*esd.fTracks)),
151   fMuonTracks(new TClonesArray(*esd.fMuonTracks)),
152   fPmdTracks(new TClonesArray(*esd.fPmdTracks)),
153   fTrdTracks(new TClonesArray(*esd.fTrdTracks)),
154   fV0s(new TClonesArray(*esd.fV0s)),  
155   fCascades(new TClonesArray(*esd.fCascades)),
156   fKinks(new TClonesArray(*esd.fKinks)),
157   fCaloClusters(new TClonesArray(*esd.fCaloClusters)),
158   fEMCALCells(new AliESDCaloCells(*esd.fEMCALCells)),
159   fPHOSCells(new AliESDCaloCells(*esd.fPHOSCells)),
160   fErrorLogs(new TClonesArray(*esd.fErrorLogs)),
161   fESDOld(new AliESD(*esd.fESDOld)),
162   fESDFriendOld(new AliESDfriend(*esd.fESDFriendOld)),
163   fConnected(esd.fConnected),
164   fUseOwnList(esd.fUseOwnList),
165   fEMCALClusters(esd.fEMCALClusters), 
166   fFirstEMCALCluster(esd.fFirstEMCALCluster),
167   fPHOSClusters(esd.fPHOSClusters), 
168   fFirstPHOSCluster(esd.fFirstPHOSCluster)
169
170 {
171   // CKB init in the constructor list and only add here ...
172   AddObject(fESDRun);
173   AddObject(fHeader);
174   AddObject(fESDZDC);
175   AddObject(fESDFMD);
176   AddObject(fESDVZERO);
177   AddObject(fESDTZERO);
178   AddObject(fTPCVertex);
179   AddObject(fSPDVertex);
180   AddObject(fPrimaryVertex);
181   AddObject(fSPDMult);
182   AddObject(fPHOSTrigger);
183   AddObject(fEMCALTrigger);
184   AddObject(fTracks);
185   AddObject(fMuonTracks);
186   AddObject(fPmdTracks);
187   AddObject(fTrdTracks);
188   AddObject(fV0s);
189   AddObject(fCascades);
190   AddObject(fKinks);
191   AddObject(fCaloClusters);
192   AddObject(fEMCALCells);
193   AddObject(fPHOSCells);
194   AddObject(fErrorLogs);
195   AddObject(fESDACORDE);
196
197   GetStdContent();
198
199 }
200
201 //______________________________________________________________________________
202 AliESDEvent & AliESDEvent::operator=(const AliESDEvent& source) {
203
204   // Assignment operator
205
206   if(&source == this) return *this;
207   AliVEvent::operator=(source);
208
209   // This assumes that the list is already created
210   // and that the virtual void Copy(Tobject&) function
211   // is correctly implemented in the derived class
212   // otherwise only TObject::Copy() will be used
213
214
215
216   if((fESDObjects->GetSize()==0)&&(source.fESDObjects->GetSize()>=kESDListN)){
217     // We cover the case that we do not yet have the 
218     // standard content but the source has it
219     CreateStdContent();
220   }
221
222   TIter next(source.GetList());
223   TObject *its = 0;
224   TString name;
225   while ((its = next())) {
226     name.Form("%s", its->GetName());
227     TObject *mine = fESDObjects->FindObject(name.Data());
228     if(!mine){
229       TClass* pClass=TClass::GetClass(its->ClassName());
230       if (!pClass) {
231         AliWarning(Form("Can not find class description for entry %s (%s)\n",
232                         its->ClassName(), name.Data()));
233         continue;
234       }
235
236       mine=(TObject*)pClass->New();
237       if(!mine){
238       // not in this: can be added to list
239         AliWarning(Form("%s:%d Could not find %s for copying \n",
240                         (char*)__FILE__,__LINE__,name.Data()));
241         continue;
242       }  
243       if(mine->InheritsFrom("TNamed")){
244         ((TNamed*)mine)->SetName(name);
245       }
246       else if(mine->InheritsFrom("TCollection")){
247         if(mine->InheritsFrom("TClonesArray"))
248           dynamic_cast<TClonesArray*>(mine)->SetClass(dynamic_cast<TClonesArray*>(its)->GetClass());
249         dynamic_cast<TCollection*>(mine)->SetName(name);
250       }
251       AliDebug(1, Form("adding object %s of type %s", mine->GetName(), mine->ClassName()));
252       AddObject(mine);
253     }  
254    
255     if(!its->InheritsFrom("TCollection")){
256       // simple objects
257       its->Copy(*mine);
258     }
259     else if(its->InheritsFrom("TClonesArray")){
260       // Create or expand the tclonesarray pointers
261       // so we can directly copy to the object
262       TClonesArray *its_tca = (TClonesArray*)its;
263       TClonesArray *mine_tca = (TClonesArray*)mine;
264
265       // this leaves the capacity of the TClonesArray the same
266       // except for a factor of 2 increase when size > capacity
267       // does not release any memory occupied by the tca
268       mine_tca->ExpandCreate(its_tca->GetEntriesFast());
269       for(int i = 0;i < its_tca->GetEntriesFast();++i){
270         // copy 
271         TObject *mine_tca_obj = mine_tca->At(i);
272         TObject *its_tca_obj = its_tca->At(i);
273         // no need to delete first
274         // pointers within the class should be handled by Copy()...
275         // Can there be Empty slots?
276         its_tca_obj->Copy(*mine_tca_obj);
277       }
278     }
279     else{
280       AliWarning(Form("%s:%d cannot copy TCollection \n",
281                       (char*)__FILE__,__LINE__));
282     }
283   }
284
285   fConnected = source.fConnected;
286   fUseOwnList = source.fUseOwnList;
287   fEMCALClusters = source.fEMCALClusters;
288   fFirstEMCALCluster = source.fFirstEMCALCluster;
289   fPHOSClusters = source.fPHOSClusters;
290   fFirstPHOSCluster = source.fFirstPHOSCluster;
291
292
293   return *this;
294
295 }
296
297
298 //______________________________________________________________________________
299 AliESDEvent::~AliESDEvent()
300 {
301   //
302   // Standard destructor
303   //
304
305   // everthing on the list gets deleted automatically
306
307   
308   if(fESDObjects&&!fConnected)
309     {
310       delete fESDObjects;
311       fESDObjects = 0;
312     }
313
314   
315 }
316
317 void AliESDEvent::Copy(TObject &obj) const {
318
319   // interface to TOBject::Copy
320   // Copies the content of this into obj!
321   // bascially obj = *this
322
323   if(this==&obj)return;
324   AliESDEvent *robj = dynamic_cast<AliESDEvent*>(&obj);
325   if(!robj)return; // not an AliESEvent
326   *robj = *this;
327   return;
328 }
329
330 //______________________________________________________________________________
331 void AliESDEvent::Reset()
332 {
333
334   // Handle the cases
335   // Std content + Non std content
336
337   // Reset the standard contents
338   ResetStdContent(); 
339
340   //  reset for the old data without AliESDEvent...
341   if(fESDOld)fESDOld->Reset();
342   if(fESDFriendOld){
343     fESDFriendOld->~AliESDfriend();
344     new (fESDFriendOld) AliESDfriend();
345   }
346   // 
347
348   if(fESDObjects->GetSize()>kESDListN){
349     // we have non std content
350     // this also covers esdfriends
351     for(int i = kESDListN;i < fESDObjects->GetSize();++i){
352       TObject *pObject = fESDObjects->At(i);
353       // TClonesArrays
354       if(pObject->InheritsFrom(TClonesArray::Class())){
355         ((TClonesArray*)pObject)->Delete();
356       }
357       else if(!pObject->InheritsFrom(TCollection::Class())){
358         ResetWithPlacementNew(pObject);
359       }
360       else{
361         AliWarning(Form("No reset for %s (%s)\n",
362                         pObject->ClassName()));
363       }
364     }
365   }
366
367 }
368
369 Bool_t AliESDEvent::ResetWithPlacementNew(TObject *pObject){
370   Long_t dtoronly = TObject::GetDtorOnly();
371   TClass *pClass = TClass::GetClass(pObject->ClassName()); 
372   TObject::SetDtorOnly(pObject);
373   delete pObject;
374   // Recreate with placement new
375   pClass->New(pObject);
376   // Restore the state.
377   TObject::SetDtorOnly((void*)dtoronly);
378   return kTRUE;
379 }
380
381 void AliESDEvent::ResetStdContent()
382 {
383   // Reset the standard contents
384   if(fESDRun) fESDRun->Reset();
385   if(fHeader) fHeader->Reset();
386   if(fESDZDC) fESDZDC->Reset();
387   if(fESDFMD) {
388     fESDFMD->Clear();
389   }
390   if(fESDVZERO){
391     // reset by callin d'to /c'tor keep the pointer
392     fESDVZERO->~AliESDVZERO();
393     new (fESDVZERO) AliESDVZERO();
394   }  
395   if(fESDACORDE){
396     fESDACORDE->~AliESDACORDE();
397     new (fESDACORDE) AliESDACORDE();    
398   } 
399   if(fESDTZERO) fESDTZERO->Reset(); 
400   // CKB no clear/reset implemented
401   if(fTPCVertex){
402     fTPCVertex->~AliESDVertex();
403     new (fTPCVertex) AliESDVertex();
404     fTPCVertex->SetName(fgkESDListName[kTPCVertex]);
405   }
406   if(fSPDVertex){
407     fSPDVertex->~AliESDVertex();
408     new (fSPDVertex) AliESDVertex();
409     fSPDVertex->SetName(fgkESDListName[kSPDVertex]);
410   }
411   if(fPrimaryVertex){
412     fPrimaryVertex->~AliESDVertex();
413     new (fPrimaryVertex) AliESDVertex();
414     fPrimaryVertex->SetName(fgkESDListName[kPrimaryVertex]);
415   }
416   if(fSPDMult){
417     fSPDMult->~AliMultiplicity();
418     new (fSPDMult) AliMultiplicity();
419   }
420   if(fPHOSTrigger)fPHOSTrigger->Reset(); 
421   if(fEMCALTrigger)fEMCALTrigger->Reset(); 
422   if(fTracks)fTracks->Delete();
423   if(fMuonTracks)fMuonTracks->Delete();
424   if(fPmdTracks)fPmdTracks->Delete();
425   if(fTrdTracks)fTrdTracks->Delete();
426   if(fV0s)fV0s->Delete();
427   if(fCascades)fCascades->Delete();
428   if(fKinks)fKinks->Delete();
429   if(fCaloClusters)fCaloClusters->Delete();
430   if(fPHOSCells)fPHOSCells->DeleteContainer();
431   if(fEMCALCells)fEMCALCells->DeleteContainer();
432   if(fErrorLogs) fErrorLogs->Delete();
433
434   // don't reset fconnected fConnected and the list
435
436   fEMCALClusters=0; 
437   fFirstEMCALCluster=-1; 
438   fPHOSClusters=0; 
439   fFirstPHOSCluster=-1; 
440 }
441
442
443 Int_t AliESDEvent::AddV0(const AliESDv0 *v) {
444   //
445   // Add V0
446   //
447   TClonesArray &fv = *fV0s;
448   Int_t idx=fV0s->GetEntriesFast();
449   new(fv[idx]) AliESDv0(*v);
450   return idx;
451 }  
452
453 //______________________________________________________________________________
454 void AliESDEvent::Print(Option_t *) const 
455 {
456   //
457   // Print header information of the event
458   //
459   printf("ESD run information\n");
460   printf("Event # in file %d Bunch crossing # %d Orbit # %d Period # %d Run # %d Trigger %lld Magnetic field %f \n",
461          GetEventNumberInFile(),
462          GetBunchCrossNumber(),
463          GetOrbitNumber(),
464          GetPeriodNumber(),
465          GetRunNumber(),
466          GetTriggerMask(),
467          GetMagneticField() );
468   printf("Vertex: (%.4f +- %.4f, %.4f +- %.4f, %.4f +- %.4f) cm\n",
469            fPrimaryVertex->GetXv(), fPrimaryVertex->GetXRes(),
470            fPrimaryVertex->GetYv(), fPrimaryVertex->GetYRes(),
471            fPrimaryVertex->GetZv(), fPrimaryVertex->GetZRes());
472     printf("Mean vertex in RUN: X=%.4f Y=%.4f cm\n",
473            GetDiamondX(),GetDiamondY());
474     printf("SPD Multiplicity. Number of tracklets %d \n",
475            fSPDMult->GetNumberOfTracklets());
476   printf("Number of tracks: \n");
477   printf("                 charged   %d\n", GetNumberOfTracks());
478   printf("                 muon      %d\n", GetNumberOfMuonTracks());
479   printf("                 pmd       %d\n", GetNumberOfPmdTracks());
480   printf("                 trd       %d\n", GetNumberOfTrdTracks());
481   printf("                 v0        %d\n", GetNumberOfV0s());
482   printf("                 cascades  %d\n", GetNumberOfCascades());
483   printf("                 kinks     %d\n", GetNumberOfKinks());
484   if(fPHOSCells)printf("                 PHOSCells %d\n", fPHOSCells->GetNumberOfCells());
485   else printf("                 PHOSCells not in the Event\n");
486   if(fEMCALCells)printf("                 EMCALCells %d\n", fEMCALCells->GetNumberOfCells());
487   else printf("                 EMCALCells not in the Event\n");
488   printf("                 CaloClusters %d\n", GetNumberOfCaloClusters());
489   printf("                 phos      %d\n", GetNumberOfPHOSClusters());
490   printf("                 emcal     %d\n", GetNumberOfEMCALClusters());
491   printf("                 FMD       %s\n", (fESDFMD ? "yes" : "no"));
492   printf("                 VZERO     %s\n", (fESDVZERO ? "yes" : "no"));
493
494   return;
495 }
496
497 void AliESDEvent::SetESDfriend(const AliESDfriend *ev) const {
498   //
499   // Attaches the complementary info to the ESD
500   //
501   if (!ev) return;
502
503   // to be sure that we set the tracks also
504   // in case of old esds 
505   // if(fESDOld)CopyFromOldESD();
506
507   Int_t ntrk=ev->GetNumberOfTracks();
508  
509   for (Int_t i=0; i<ntrk; i++) {
510     const AliESDfriendTrack *f=ev->GetTrack(i);
511     GetTrack(i)->SetFriendTrack(f);
512   }
513 }
514
515 Bool_t  AliESDEvent::RemoveKink(Int_t rm) const {
516   // ---------------------------------------------------------
517   // Remove a kink candidate and references to it from ESD,
518   // if this candidate does not come from a reconstructed decay
519   // Not yet implemented...
520   // ---------------------------------------------------------
521   Int_t last=GetNumberOfKinks()-1;
522   if ((rm<0)||(rm>last)) return kFALSE;
523
524   return kTRUE;
525 }
526
527 Bool_t  AliESDEvent::RemoveV0(Int_t rm) const {
528   // ---------------------------------------------------------
529   // Remove a V0 candidate and references to it from ESD,
530   // if this candidate does not come from a reconstructed decay
531   // ---------------------------------------------------------
532   Int_t last=GetNumberOfV0s()-1;
533   if ((rm<0)||(rm>last)) return kFALSE;
534
535   AliESDv0 *v0=GetV0(rm);
536   Int_t idxP=v0->GetPindex(), idxN=v0->GetNindex();
537
538   v0=GetV0(last);
539   Int_t lastIdxP=v0->GetPindex(), lastIdxN=v0->GetNindex();
540
541   Int_t used=0;
542
543   // Check if this V0 comes from a reconstructed decay
544   Int_t ncs=GetNumberOfCascades();
545   for (Int_t n=0; n<ncs; n++) {
546     AliESDcascade *cs=GetCascade(n);
547
548     Int_t csIdxP=cs->GetPindex();
549     Int_t csIdxN=cs->GetNindex();
550
551     if (idxP==csIdxP)
552        if (idxN==csIdxN) return kFALSE;
553
554     if (csIdxP==lastIdxP)
555        if (csIdxN==lastIdxN) used++;
556   }
557
558   //Replace the removed V0 with the last V0 
559   TClonesArray &a=*fV0s;
560   delete a.RemoveAt(rm);
561
562   if (rm==last) return kTRUE;
563
564   //v0 is pointing to the last V0 candidate... 
565   new (a[rm]) AliESDv0(*v0);
566   delete a.RemoveAt(last);
567
568   if (!used) return kTRUE;
569   
570
571   // Remap the indices of the daughters of reconstructed decays
572   for (Int_t n=0; n<ncs; n++) {
573     AliESDcascade *cs=GetCascade(n);
574
575
576     Int_t csIdxP=cs->GetPindex();
577     Int_t csIdxN=cs->GetNindex();
578
579     if (csIdxP==lastIdxP)
580       if (csIdxN==lastIdxN) {
581          cs->AliESDv0::SetIndex(1,idxP);
582          cs->AliESDv0::SetIndex(0,idxN);
583          used--;
584          if (!used) return kTRUE;
585       }
586   }
587
588   return kTRUE;
589 }
590
591 Bool_t  AliESDEvent::RemoveTrack(Int_t rm) const {
592   // ---------------------------------------------------------
593   // Remove a track and references to it from ESD,
594   // if this track does not come from a reconstructed decay
595   // ---------------------------------------------------------
596   Int_t last=GetNumberOfTracks()-1;
597   if ((rm<0)||(rm>last)) return kFALSE;
598
599   Int_t used=0;
600
601   // Check if this track comes from the reconstructed primary vertices
602   if (fTPCVertex && fTPCVertex->GetStatus()) {
603      UShort_t *primIdx=fTPCVertex->GetIndices();
604      Int_t n=fTPCVertex->GetNIndices();
605      while (n--) {
606        Int_t idx=Int_t(primIdx[n]);
607        if (rm==idx) return kFALSE;
608        if (idx==last) used++; 
609      }
610   }
611   if (fPrimaryVertex && fPrimaryVertex->GetStatus()) {
612      UShort_t *primIdx=fPrimaryVertex->GetIndices();
613      Int_t n=fPrimaryVertex->GetNIndices();
614      while (n--) {
615        Int_t idx=Int_t(primIdx[n]);
616        if (rm==idx) return kFALSE;
617        if (idx==last) used++; 
618      }
619   }
620   
621   // Check if this track comes from a reconstructed decay
622   Int_t nv0=GetNumberOfV0s();
623   for (Int_t n=0; n<nv0; n++) {
624     AliESDv0 *v0=GetV0(n);
625
626     Int_t idx=v0->GetNindex();
627     if (rm==idx) return kFALSE;
628     if (idx==last) used++;
629
630     idx=v0->GetPindex();
631     if (rm==idx) return kFALSE;
632     if (idx==last) used++;
633   }
634
635   Int_t ncs=GetNumberOfCascades();
636   for (Int_t n=0; n<ncs; n++) {
637     AliESDcascade *cs=GetCascade(n);
638
639     Int_t idx=cs->GetIndex();
640     if (rm==idx) return kFALSE;
641     if (idx==last) used++;
642
643     AliESDv0 *v0=cs;
644     idx=v0->GetNindex();
645     if (rm==idx) return kFALSE;
646     if (idx==last) used++;
647
648     idx=v0->GetPindex();
649     if (rm==idx) return kFALSE;
650     if (idx==last) used++;
651   }
652
653   Int_t nkn=GetNumberOfKinks();
654   for (Int_t n=0; n<nkn; n++) {
655     AliESDkink *kn=GetKink(n);
656
657     Int_t idx=kn->GetIndex(0);
658     if (rm==idx) return kFALSE;
659     if (idx==last) used++;
660
661     idx=kn->GetIndex(1);
662     if (rm==idx) return kFALSE;
663     if (idx==last) used++;
664   }
665
666
667   //Replace the removed track with the last track 
668   TClonesArray &a=*fTracks;
669   delete a.RemoveAt(rm);
670
671   if (rm==last) return kTRUE;
672
673   AliESDtrack *t=GetTrack(last);
674   t->SetID(rm);
675   new (a[rm]) AliESDtrack(*t);
676   delete a.RemoveAt(last);
677
678
679   if (!used) return kTRUE;
680   
681
682   // Remap the indices of the tracks used for the primary vertex reconstruction
683   if (fTPCVertex && fTPCVertex->GetStatus()) {
684      UShort_t *primIdx=fTPCVertex->GetIndices();
685      Int_t n=fTPCVertex->GetNIndices();
686      while (n--) {
687        Int_t idx=Int_t(primIdx[n]);
688        if (idx==last) {
689           primIdx[n]=Short_t(rm); 
690           used--;
691           if (!used) return kTRUE;
692        }
693      }
694   }  
695   if (fPrimaryVertex && fPrimaryVertex->GetStatus()) {
696      UShort_t *primIdx=fPrimaryVertex->GetIndices();
697      Int_t n=fPrimaryVertex->GetNIndices();
698      while (n--) {
699        Int_t idx=Int_t(primIdx[n]);
700        if (idx==last) {
701           primIdx[n]=Short_t(rm); 
702           used--;
703           if (!used) return kTRUE;
704        }
705      }
706   }  
707
708   // Remap the indices of the daughters of reconstructed decays
709   for (Int_t n=0; n<nv0; n++) {
710     AliESDv0 *v0=GetV0(n);
711     if (v0->GetIndex(0)==last) {
712        v0->SetIndex(0,rm);
713        used--;
714        if (!used) return kTRUE;
715     }
716     if (v0->GetIndex(1)==last) {
717        v0->SetIndex(1,rm);
718        used--;
719        if (!used) return kTRUE;
720     }
721   }
722
723   for (Int_t n=0; n<ncs; n++) {
724     AliESDcascade *cs=GetCascade(n);
725     if (cs->GetIndex()==last) {
726        cs->SetIndex(rm);
727        used--;
728        if (!used) return kTRUE;
729     }
730     AliESDv0 *v0=cs;
731     if (v0->GetIndex(0)==last) {
732        v0->SetIndex(0,rm);
733        used--;
734        if (!used) return kTRUE;
735     }
736     if (v0->GetIndex(1)==last) {
737        v0->SetIndex(1,rm);
738        used--;
739        if (!used) return kTRUE;
740     }
741   }
742
743   for (Int_t n=0; n<nkn; n++) {
744     AliESDkink *kn=GetKink(n);
745     if (kn->GetIndex(0)==last) {
746        kn->SetIndex(rm,0);
747        used--;
748        if (!used) return kTRUE;
749     }
750     if (kn->GetIndex(1)==last) {
751        kn->SetIndex(rm,1);
752        used--;
753        if (!used) return kTRUE;
754     }
755   }
756
757   return kTRUE;
758 }
759
760
761 Bool_t AliESDEvent::Clean(Float_t *cleanPars) {
762   //
763   // Remove the data which are not needed for the physics analysis.
764   //
765   // 1) Cleaning the V0 candidates
766   //    ---------------------------
767   //    If the cosine of the V0 pointing angle "csp" and 
768   //    the DCA between the daughter tracks "dca" does not satisfy 
769   //    the conditions 
770   //
771   //     csp > cleanPars[1] + dca/cleanPars[0]*(1.- cleanPars[1])
772   //
773   //    an attempt to remove this V0 candidate from ESD is made.
774   //
775   //    The V0 candidate gets removed if it does not belong to any 
776   //    recosntructed cascade decay
777   //
778   //    12.11.2007, optimal values: cleanPars[0]=0.5, cleanPars[1]=0.999
779   //
780   // 2) Cleaning the tracks
781   //    ----------------------
782   //    If track's transverse parameter is larger than cleanPars[2]
783   //                       OR
784   //    track's longitudinal parameter is larger than cleanPars[3]
785   //    an attempt to remove this track from ESD is made.
786   //
787   //    The track gets removed if it does not come 
788   //    from a reconstructed decay
789   //
790   Bool_t rc=kFALSE;
791
792   Float_t dcaMax=cleanPars[0];
793   Float_t cspMin=cleanPars[1];
794
795   Int_t nV0s=GetNumberOfV0s();
796   for (Int_t i=nV0s-1; i>=0; i--) {
797     AliESDv0 *v0=GetV0(i);
798
799     Float_t dca=v0->GetDcaV0Daughters();
800     Float_t csp=v0->GetV0CosineOfPointingAngle();
801     Float_t cspcut=cspMin + dca/dcaMax*(1.-cspMin);
802     if (csp > cspcut) continue;
803     if (RemoveV0(i)) rc=kTRUE;
804   }
805
806
807   Float_t dmax=cleanPars[2], zmax=cleanPars[3];
808
809   const AliESDVertex *vertex=GetPrimaryVertexSPD();
810   Bool_t vtxOK=vertex->GetStatus();
811   
812   Int_t nTracks=GetNumberOfTracks();
813   for (Int_t i=nTracks-1; i>=0; i--) {
814     AliESDtrack *track=GetTrack(i);
815     Float_t xy,z; track->GetImpactParameters(xy,z);
816     if ((TMath::Abs(xy) > dmax) || (vtxOK && (TMath::Abs(z) > zmax))) {
817       if (RemoveTrack(i)) rc=kTRUE;
818     }
819   }
820
821   return rc;
822 }
823
824 Int_t  AliESDEvent::AddTrack(const AliESDtrack *t) 
825 {
826     // Add track
827     TClonesArray &ftr = *fTracks;
828     AliESDtrack * track = new(ftr[fTracks->GetEntriesFast()])AliESDtrack(*t);
829     track->SetID(fTracks->GetEntriesFast()-1);
830     return  track->GetID();    
831 }
832
833  void AliESDEvent::AddMuonTrack(const AliESDMuonTrack *t) 
834 {
835     TClonesArray &fmu = *fMuonTracks;
836     new(fmu[fMuonTracks->GetEntriesFast()]) AliESDMuonTrack(*t);
837 }
838
839 void AliESDEvent::AddPmdTrack(const AliESDPmdTrack *t) 
840 {
841   TClonesArray &fpmd = *fPmdTracks;
842   new(fpmd[fPmdTracks->GetEntriesFast()]) AliESDPmdTrack(*t);
843 }
844
845 void AliESDEvent::AddTrdTrack(const AliESDTrdTrack *t) 
846 {
847   TClonesArray &ftrd = *fTrdTracks;
848   new(ftrd[fTrdTracks->GetEntriesFast()]) AliESDTrdTrack(*t);
849 }
850
851
852
853
854 Int_t AliESDEvent::AddKink(const AliESDkink *c) 
855 {
856     // Add kink
857     TClonesArray &fk = *fKinks;
858     AliESDkink * kink = new(fk[fKinks->GetEntriesFast()]) AliESDkink(*c);
859     kink->SetID(fKinks->GetEntriesFast()); // CKB different from the other imps..
860     return fKinks->GetEntriesFast()-1;
861 }
862
863
864 void AliESDEvent::AddCascade(const AliESDcascade *c) 
865 {
866   TClonesArray &fc = *fCascades;
867   new(fc[fCascades->GetEntriesFast()]) AliESDcascade(*c);
868 }
869
870
871 Int_t AliESDEvent::AddCaloCluster(const AliESDCaloCluster *c) 
872 {
873     // Add calocluster
874     TClonesArray &fc = *fCaloClusters;
875     AliESDCaloCluster *clus = new(fc[fCaloClusters->GetEntriesFast()]) AliESDCaloCluster(*c);
876     clus->SetID(fCaloClusters->GetEntriesFast()-1);
877     return fCaloClusters->GetEntriesFast()-1;
878   }
879
880
881 void  AliESDEvent::AddRawDataErrorLog(const AliRawDataErrorLog *log) const {
882   TClonesArray &errlogs = *fErrorLogs;
883   new(errlogs[errlogs.GetEntriesFast()])  AliRawDataErrorLog(*log);
884 }
885
886 void  AliESDEvent::SetPrimaryVertexTPC(const AliESDVertex *vertex) 
887 {
888   // Set the TPC vertex
889   // use already allocated space
890   if(fTPCVertex){
891     *fTPCVertex = *vertex;
892     fTPCVertex->SetName(fgkESDListName[kTPCVertex]);
893   }
894 }
895
896 void  AliESDEvent::SetPrimaryVertexSPD(const AliESDVertex *vertex) 
897 {
898   // Set the SPD vertex
899   // use already allocated space
900   if(fSPDVertex){
901     *fSPDVertex = *vertex;
902     fSPDVertex->SetName(fgkESDListName[kSPDVertex]);
903   }
904 }
905
906 void  AliESDEvent::SetPrimaryVertexTracks(const AliESDVertex *vertex) 
907 {
908   // Set the primary vertex reconstructed using he ESD tracks.
909   // use already allocated space
910   if(fPrimaryVertex){
911     *fPrimaryVertex = *vertex;
912     fPrimaryVertex->SetName(fgkESDListName[kPrimaryVertex]);
913   }
914 }
915
916 const AliESDVertex * AliESDEvent::GetPrimaryVertex() const 
917 {
918   //
919   // Get the "best" available reconstructed primary vertex.
920   //
921   if(fPrimaryVertex){
922     if (fPrimaryVertex->GetStatus()) return fPrimaryVertex;
923   }
924   if(fSPDVertex){
925     if (fSPDVertex->GetStatus()) return fSPDVertex;
926   }
927   if(fTPCVertex) return fTPCVertex;
928   
929   AliWarning("No primary vertex available. Returning the \"default\"...");
930   return fSPDVertex;
931 }
932
933 void AliESDEvent::SetMultiplicity(const AliMultiplicity *mul) 
934 {
935   // Set the SPD Multiplicity
936   if(fSPDMult){
937     *fSPDMult = *mul;
938   }
939 }
940
941
942 void AliESDEvent::SetFMDData(AliESDFMD * obj) 
943
944   // use already allocated space
945   if(fESDFMD){
946     *fESDFMD = *obj;
947   }
948 }
949
950 void AliESDEvent::SetVZEROData(AliESDVZERO * obj)
951
952   // use already allocated space
953   if(fESDVZERO)
954     *fESDVZERO = *obj;
955 }
956
957 void AliESDEvent::SetACORDEData(AliESDACORDE * obj)
958 {
959   if(fESDACORDE)
960     *fESDACORDE = *obj;
961 }
962
963
964 void AliESDEvent::GetESDfriend(AliESDfriend *ev) const 
965 {
966   //
967   // Extracts the complementary info from the ESD
968   //
969   if (!ev) return;
970
971   Int_t ntrk=GetNumberOfTracks();
972
973   for (Int_t i=0; i<ntrk; i++) {
974     AliESDtrack *t=GetTrack(i);
975     const AliESDfriendTrack *f=t->GetFriendTrack();
976     ev->AddTrack(f);
977
978     t->ReleaseESDfriendTrack();// Not to have two copies of "friendTrack"
979
980   }
981
982   AliESDfriend *fr = (AliESDfriend*)(const_cast<AliESDEvent*>(this)->FindListObject("AliESDfriend"));
983   if (fr) ev->SetVZEROfriend(fr->GetVZEROfriend());
984 }
985
986 void AliESDEvent::AddObject(TObject* obj) 
987 {
988   // Add an object to the list of object.
989   // Please be aware that in order to increase performance you should
990   // refrain from using TObjArrays (if possible). Use TClonesArrays, instead.
991   fESDObjects->SetOwner(kTRUE);
992   fESDObjects->AddLast(obj);
993 }
994
995
996 void AliESDEvent::GetStdContent() 
997 {
998   // set pointers for standard content
999   // get by name much safer and not a big overhead since not called very often
1000  
1001   fESDRun = (AliESDRun*)fESDObjects->FindObject(fgkESDListName[kESDRun]);
1002   fHeader = (AliESDHeader*)fESDObjects->FindObject(fgkESDListName[kHeader]);
1003   fESDZDC = (AliESDZDC*)fESDObjects->FindObject(fgkESDListName[kESDZDC]);
1004   fESDFMD = (AliESDFMD*)fESDObjects->FindObject(fgkESDListName[kESDFMD]);
1005   fESDVZERO = (AliESDVZERO*)fESDObjects->FindObject(fgkESDListName[kESDVZERO]);
1006   fESDTZERO = (AliESDTZERO*)fESDObjects->FindObject(fgkESDListName[kESDTZERO]);
1007   fTPCVertex = (AliESDVertex*)fESDObjects->FindObject(fgkESDListName[kTPCVertex]);
1008   fSPDVertex = (AliESDVertex*)fESDObjects->FindObject(fgkESDListName[kSPDVertex]);
1009   fPrimaryVertex = (AliESDVertex*)fESDObjects->FindObject(fgkESDListName[kPrimaryVertex]);
1010   fSPDMult =       (AliMultiplicity*)fESDObjects->FindObject(fgkESDListName[kSPDMult]);
1011   fPHOSTrigger = (AliESDCaloTrigger*)fESDObjects->FindObject(fgkESDListName[kPHOSTrigger]);
1012   fEMCALTrigger = (AliESDCaloTrigger*)fESDObjects->FindObject(fgkESDListName[kEMCALTrigger]);
1013   fTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kTracks]);
1014   fMuonTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kMuonTracks]);
1015   fPmdTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kPmdTracks]);
1016   fTrdTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kTrdTracks]);
1017   fV0s = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kV0s]);
1018   fCascades = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kCascades]);
1019   fKinks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kKinks]);
1020   fCaloClusters = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kCaloClusters]);
1021   fEMCALCells = (AliESDCaloCells*)fESDObjects->FindObject(fgkESDListName[kEMCALCells]);
1022   fPHOSCells = (AliESDCaloCells*)fESDObjects->FindObject(fgkESDListName[kPHOSCells]);
1023   fErrorLogs = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kErrorLogs]);
1024   fESDACORDE = (AliESDACORDE*)fESDObjects->FindObject(fgkESDListName[kESDACORDE]);
1025
1026 }
1027
1028 void AliESDEvent::SetStdNames(){
1029   // Set the names of the standard contents
1030   // 
1031   if(fESDObjects->GetEntries()>=kESDListN){
1032     for(int i = 0;i < fESDObjects->GetEntries() && i<kESDListN;i++){
1033       TObject *fObj = fESDObjects->At(i);
1034       if(fObj->InheritsFrom("TNamed")){
1035         ((TNamed*)fObj)->SetName(fgkESDListName[i]);
1036       }
1037       else if(fObj->InheritsFrom("TClonesArray")){
1038         ((TClonesArray*)fObj)->SetName(fgkESDListName[i]);
1039       }
1040     }
1041   }
1042   else{
1043      AliWarning("Std Entries missing");
1044   }
1045
1046
1047
1048 void AliESDEvent::CreateStdContent(Bool_t bUseThisList){
1049   fUseOwnList = bUseThisList;
1050   CreateStdContent();
1051 }
1052
1053 void AliESDEvent::CreateStdContent() 
1054 {
1055   // create the standard AOD content and set pointers
1056
1057   // create standard objects and add them to the TList of objects
1058   AddObject(new AliESDRun());
1059   AddObject(new AliESDHeader());
1060   AddObject(new AliESDZDC());
1061   AddObject(new AliESDFMD());
1062   AddObject(new AliESDVZERO());
1063   AddObject(new AliESDTZERO());
1064   AddObject(new AliESDVertex());
1065   AddObject(new AliESDVertex());
1066   AddObject(new AliESDVertex());
1067   AddObject(new AliMultiplicity());
1068   AddObject(new AliESDCaloTrigger());
1069   AddObject(new AliESDCaloTrigger());
1070   AddObject(new TClonesArray("AliESDtrack",0));
1071   AddObject(new TClonesArray("AliESDMuonTrack",0));
1072   AddObject(new TClonesArray("AliESDPmdTrack",0));
1073   AddObject(new TClonesArray("AliESDTrdTrack",0));
1074   AddObject(new TClonesArray("AliESDv0",0));
1075   AddObject(new TClonesArray("AliESDcascade",0));
1076   AddObject(new TClonesArray("AliESDkink",0));
1077   AddObject(new TClonesArray("AliESDCaloCluster",0));
1078   AddObject(new AliESDCaloCells());
1079   AddObject(new AliESDCaloCells());
1080   AddObject(new TClonesArray("AliRawDataErrorLog",0));
1081   AddObject(new AliESDACORDE()); 
1082
1083   // check the order of the indices against enum...
1084
1085   // set names
1086   SetStdNames();
1087   // read back pointers
1088   GetStdContent();
1089 }
1090
1091 TObject* AliESDEvent::FindListObject(const char *name){
1092 //
1093 // Find object with name "name" in the list of branches
1094 //
1095   if(fESDObjects){
1096     return fESDObjects->FindObject(name);
1097   }
1098   return 0;
1099
1100
1101 Int_t AliESDEvent::GetPHOSClusters(TRefArray *clusters) const
1102 {
1103   // fills the provided TRefArray with all found phos clusters
1104   
1105   clusters->Clear();
1106   
1107   AliESDCaloCluster *cl = 0;
1108   for (Int_t i = 0; i < GetNumberOfCaloClusters(); i++) {
1109     
1110     if ( (cl = GetCaloCluster(i)) ) {
1111       if (cl->IsPHOS()){
1112         clusters->Add(cl);
1113         AliDebug(1,Form("IsPHOS cluster %d Size: %d \n",i,clusters->GetEntriesFast()));
1114       }
1115     }
1116   }
1117   return clusters->GetEntriesFast();
1118 }
1119
1120 Int_t AliESDEvent::GetEMCALClusters(TRefArray *clusters) const
1121 {
1122   // fills the provided TRefArray with all found emcal clusters
1123
1124   clusters->Clear();
1125
1126   AliESDCaloCluster *cl = 0;
1127   for (Int_t i = 0; i < GetNumberOfCaloClusters(); i++) {
1128
1129     if ( (cl = GetCaloCluster(i)) ) {
1130       if (cl->IsEMCAL()){
1131         clusters->Add(cl);
1132         AliDebug(1,Form("IsEMCAL cluster %d Size: %d \n",i,clusters->GetEntriesFast()));
1133       }
1134     }
1135   }
1136   return clusters->GetEntriesFast();
1137 }
1138
1139 void AliESDEvent::WriteToTree(TTree* tree) const {
1140   // Book the branches as in TTree::Branch(TCollection*)
1141   // but add a "." at the end of top level branches which are
1142   // not a TClonesArray
1143
1144
1145   TString branchname;
1146   TIter next(fESDObjects);
1147   const Int_t kSplitlevel = 99; // default value in TTree::Branch()
1148   const Int_t kBufsize = 32000; // default value in TTree::Branch()
1149   TObject *obj = 0;
1150
1151   while ((obj = next())) {
1152     branchname.Form("%s", obj->GetName());
1153     if(branchname.CompareTo("AliESDfriend")==0)branchname = "ESDfriend.";
1154     if ((kSplitlevel > 1) &&  !obj->InheritsFrom(TClonesArray::Class())) {
1155       if(!branchname.EndsWith("."))branchname += ".";
1156     }
1157     if (!tree->FindBranch(branchname)) {
1158       tree->Bronch(branchname, obj->ClassName(), fESDObjects->GetObjectRef(obj),
1159                    kBufsize, kSplitlevel - 1);
1160     }
1161   }
1162 }
1163
1164
1165 void AliESDEvent::ReadFromTree(TTree *tree, Option_t* /*opt*/){
1166 //
1167 // Connect the ESDEvent to a tree
1168 //
1169   if(!tree){
1170     AliWarning("AliESDEvent::ReadFromTree() Zero Pointer to Tree \n");
1171     return;
1172   }
1173   // load the TTree
1174   if(!tree->GetTree())tree->LoadTree(0);
1175
1176   // if we find the "ESD" branch on the tree we do have the old structure
1177   if(tree->GetBranch("ESD")) {
1178     char ** address  = (char **)(tree->GetBranch("ESD")->GetAddress());
1179     // do we have the friend branch
1180     TBranch * esdFB = tree->GetBranch("ESDfriend.");
1181     char ** addressF = 0;
1182     if(esdFB)addressF = (char **)(esdFB->GetAddress());
1183     if (!address) {
1184       AliInfo("AliESDEvent::ReadFromTree() Reading old Tree");
1185       tree->SetBranchAddress("ESD",       &fESDOld);
1186       if(esdFB){
1187         tree->SetBranchAddress("ESDfriend.",&fESDFriendOld);
1188       }
1189     } else {
1190       AliInfo("AliESDEvent::ReadFromTree() Reading old Tree");
1191       AliInfo("Branch already connected. Using existing branch address.");
1192       fESDOld       = (AliESD*)       (*address);
1193       // addressF can still be 0, since branch needs to switched on
1194       if(addressF)fESDFriendOld = (AliESDfriend*) (*addressF);
1195     }
1196                                        
1197     //  have already connected the old ESD structure... ?
1198     // reuse also the pointer of the AlliESDEvent
1199     // otherwise create new ones
1200     TList* connectedList = (TList*) (tree->GetUserInfo()->FindObject("ESDObjectsConnectedToTree"));
1201   
1202     if(connectedList){
1203       // If connected use the connected list of objects
1204       if(fESDObjects!= connectedList){
1205         // protect when called twice 
1206         fESDObjects->Delete();
1207         fESDObjects = connectedList;
1208       }
1209       GetStdContent(); 
1210
1211       
1212       // The pointer to the friend changes when called twice via InitIO
1213       // since AliESDEvent is deleted
1214       TObject* oldf = FindListObject("AliESDfriend");
1215       TObject* newf = 0;
1216       if(addressF){
1217         newf = (TObject*)*addressF;
1218       }
1219       if(newf!=0&&oldf!=newf){
1220         // remove the old reference
1221         // Should we also delete it? Or is this handled in TTree I/O
1222         // since it is created by the first SetBranchAddress
1223         fESDObjects->Remove(oldf);
1224         // add the new one 
1225         fESDObjects->Add(newf);
1226       }
1227       
1228       fConnected = true;
1229       return;
1230     }
1231     // else...    
1232     CreateStdContent(); // create for copy
1233     // if we have the esdfriend add it, so we always can access it via the userinfo
1234     if(fESDFriendOld)AddObject(fESDFriendOld);
1235     // we are not owner of the list objects 
1236     // must not delete it
1237     fESDObjects->SetOwner(kFALSE);
1238     fESDObjects->SetName("ESDObjectsConnectedToTree");
1239     tree->GetUserInfo()->Add(fESDObjects);
1240     fConnected = true;
1241     return;
1242   }
1243   
1244
1245     delete fESDOld;
1246     fESDOld = 0;
1247   // Try to find AliESDEvent
1248   AliESDEvent *esdEvent = 0;
1249   esdEvent = (AliESDEvent*)tree->GetTree()->GetUserInfo()->FindObject("AliESDEvent");
1250   if(esdEvent){   
1251       // Check if already connected to tree
1252     TList* connectedList = (TList*) (tree->GetUserInfo()->FindObject("ESDObjectsConnectedToTree"));
1253     if (connectedList) {
1254       // If connected use the connected list if objects
1255       fESDObjects->Delete();
1256       fESDObjects = connectedList;
1257       GetStdContent(); 
1258       fConnected = true;
1259       return;
1260     }
1261
1262     // Connect to tree
1263     // prevent a memory leak when reading back the TList
1264
1265     if(!fUseOwnList){
1266       delete fESDObjects;
1267       fESDObjects = 0;
1268       // create a new TList from the UserInfo TList... 
1269       // copy constructor does not work...
1270       fESDObjects = (TList*)(esdEvent->GetList()->Clone());
1271       fESDObjects->SetOwner(kFALSE);
1272     }
1273     else if ( fESDObjects->GetEntries()==0){
1274       // at least create the std content if we want to read to our list
1275       CreateStdContent(); 
1276     }
1277
1278     // in principle
1279     // we only need new things in the list if we do no already have it..
1280     // TODO just add new entries
1281
1282     if(fESDObjects->GetEntries()<kESDListN){
1283       AliWarning(Form("AliESDEvent::ReadFromTree() TList contains less than the standard contents %d < %d \n",
1284                       fESDObjects->GetEntries(),kESDListN));
1285     }
1286     // set the branch addresses
1287     TIter next(fESDObjects);
1288     TNamed *el;
1289     while((el=(TNamed*)next())){
1290       TString bname(el->GetName());
1291       if(bname.CompareTo("AliESDfriend")==0)
1292         {
1293           // AliESDfriend does not have a name ...
1294           tree->SetBranchAddress("ESDfriend.",fESDObjects->GetObjectRef(el));
1295         }
1296       else{
1297         // check if branch exists under this Name
1298         TBranch *br = tree->GetBranch(bname.Data());
1299         if(br){
1300           tree->SetBranchAddress(bname.Data(),fESDObjects->GetObjectRef(el));
1301         }
1302         else{
1303           br = tree->GetBranch(Form("%s.",bname.Data()));
1304           if(br){
1305             tree->SetBranchAddress(Form("%s.",bname.Data()),fESDObjects->GetObjectRef(el));
1306           }
1307           else{
1308             AliWarning(Form("AliESDEvent::ReadFromTree() No Branch found with Name %s or %s.",bname.Data(),bname.Data()));
1309           }
1310
1311         }
1312       }
1313     }
1314     GetStdContent();
1315     // when reading back we are not owner of the list 
1316     // must not delete it
1317     fESDObjects->SetOwner(kFALSE);
1318     fESDObjects->SetName("ESDObjectsConnectedToTree");
1319     // we are not owner of the list objects 
1320     // must not delete it
1321     tree->GetUserInfo()->Add(fESDObjects);
1322     fConnected = true;
1323   }// no esdEvent
1324   else {
1325     // we can't get the list from the user data, create standard content
1326     // and set it by hand (no ESDfriend at the moment
1327     CreateStdContent();
1328     TIter next(fESDObjects);
1329     TNamed *el;
1330     while((el=(TNamed*)next())){
1331       TString bname(el->GetName());    
1332       TBranch *br = tree->GetBranch(bname.Data());
1333       if(br){
1334         tree->SetBranchAddress(bname.Data(),fESDObjects->GetObjectRef(el));
1335       }
1336       else{
1337         br = tree->GetBranch(Form("%s.",bname.Data()));
1338         if(br){
1339           tree->SetBranchAddress(Form("%s.",bname.Data()),fESDObjects->GetObjectRef(el));
1340         }
1341       }
1342     }
1343     GetStdContent();
1344     // when reading back we are not owner of the list 
1345     // must not delete it
1346     fESDObjects->SetOwner(kFALSE);
1347   }
1348 }
1349
1350
1351 void AliESDEvent::CopyFromOldESD()
1352 {
1353   // Method which copies over everthing from the old esd structure to the 
1354   // new  
1355   if(fESDOld){
1356     ResetStdContent();
1357      // Run
1358     SetRunNumber(fESDOld->GetRunNumber());
1359     SetPeriodNumber(fESDOld->GetPeriodNumber());
1360     SetMagneticField(fESDOld->GetMagneticField());
1361   
1362     // leave out diamond ...
1363     // SetDiamond(const AliESDVertex *vertex) { fESDRun->SetDiamond(vertex);}
1364
1365     // header
1366     SetTriggerMask(fESDOld->GetTriggerMask());
1367     SetOrbitNumber(fESDOld->GetOrbitNumber());
1368     SetTimeStamp(fESDOld->GetTimeStamp());
1369     SetEventType(fESDOld->GetEventType());
1370     SetEventNumberInFile(fESDOld->GetEventNumberInFile());
1371     SetBunchCrossNumber(fESDOld->GetBunchCrossNumber());
1372     SetTriggerCluster(fESDOld->GetTriggerCluster());
1373
1374     // ZDC
1375
1376     SetZDC(fESDOld->GetZDCN1Energy(),
1377            fESDOld->GetZDCP1Energy(),
1378            fESDOld->GetZDCEMEnergy(),
1379            0,
1380            fESDOld->GetZDCN2Energy(),
1381            fESDOld->GetZDCP2Energy(),
1382            fESDOld->GetZDCParticipants(),
1383            0);
1384
1385     // FMD
1386     
1387     if(fESDOld->GetFMDData())SetFMDData(fESDOld->GetFMDData());
1388
1389     // T0
1390
1391     SetT0zVertex(fESDOld->GetT0zVertex());
1392     SetT0(fESDOld->GetT0());
1393     //  leave amps out
1394
1395     // VZERO
1396     if (fESDOld->GetVZEROData()) SetVZEROData(fESDOld->GetVZEROData());
1397
1398     if(fESDOld->GetVertex())SetPrimaryVertexSPD(fESDOld->GetVertex());
1399
1400     if(fESDOld->GetPrimaryVertex())SetPrimaryVertexTracks(fESDOld->GetPrimaryVertex());
1401
1402     if(fESDOld->GetMultiplicity())SetMultiplicity(fESDOld->GetMultiplicity());
1403
1404     for(int i = 0;i<fESDOld->GetNumberOfTracks();i++){
1405       AddTrack(fESDOld->GetTrack(i));
1406     }
1407
1408     for(int i = 0;i<fESDOld->GetNumberOfMuonTracks();i++){
1409       AddMuonTrack(fESDOld->GetMuonTrack(i));
1410     }
1411
1412     for(int i = 0;i<fESDOld->GetNumberOfPmdTracks();i++){
1413       AddPmdTrack(fESDOld->GetPmdTrack(i));
1414     }
1415
1416     for(int i = 0;i<fESDOld->GetNumberOfTrdTracks();i++){
1417       AddTrdTrack(fESDOld->GetTrdTrack(i));
1418     }
1419
1420     for(int i = 0;i<fESDOld->GetNumberOfV0s();i++){
1421       AddV0(fESDOld->GetV0(i));
1422     }
1423
1424     for(int i = 0;i<fESDOld->GetNumberOfCascades();i++){
1425       AddCascade(fESDOld->GetCascade(i));
1426     }
1427
1428     for(int i = 0;i<fESDOld->GetNumberOfKinks();i++){
1429       AddKink(fESDOld->GetKink(i));
1430     }
1431
1432
1433     for(int i = 0;i<fESDOld->GetNumberOfCaloClusters();i++){
1434       AddCaloCluster(fESDOld->GetCaloCluster(i));
1435     }
1436
1437   }// if fesdold
1438 }
1439
1440
1441