]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliESDEvent.cxx
fNEvents initialized in the member initialization list
[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   // Check if this track is associated with a CaloCluster
667   Int_t ncl=GetNumberOfCaloClusters();
668   for (Int_t n=0; n<ncl; n++) {
669     AliESDCaloCluster *cluster=GetCaloCluster(n);
670     TArrayI *arr=cluster->GetTracksMatched();
671     Int_t s=arr->GetSize();
672     while (s--) {
673       Int_t idx=arr->At(s);
674       if (rm==idx) return kFALSE;
675       if (idx==last) used++;     
676     }
677   }
678
679
680
681   //Replace the removed track with the last track 
682   TClonesArray &a=*fTracks;
683   delete a.RemoveAt(rm);
684
685   if (rm==last) return kTRUE;
686
687   AliESDtrack *t=GetTrack(last);
688   t->SetID(rm);
689   new (a[rm]) AliESDtrack(*t);
690   delete a.RemoveAt(last);
691
692
693   if (!used) return kTRUE;
694   
695
696   // Remap the indices of the tracks used for the primary vertex reconstruction
697   if (fTPCVertex && fTPCVertex->GetStatus()) {
698      UShort_t *primIdx=fTPCVertex->GetIndices();
699      Int_t n=fTPCVertex->GetNIndices();
700      while (n--) {
701        Int_t idx=Int_t(primIdx[n]);
702        if (idx==last) {
703           primIdx[n]=Short_t(rm); 
704           used--;
705           if (!used) return kTRUE;
706        }
707      }
708   }  
709   if (fPrimaryVertex && fPrimaryVertex->GetStatus()) {
710      UShort_t *primIdx=fPrimaryVertex->GetIndices();
711      Int_t n=fPrimaryVertex->GetNIndices();
712      while (n--) {
713        Int_t idx=Int_t(primIdx[n]);
714        if (idx==last) {
715           primIdx[n]=Short_t(rm); 
716           used--;
717           if (!used) return kTRUE;
718        }
719      }
720   }  
721
722   // Remap the indices of the daughters of reconstructed decays
723   for (Int_t n=0; n<nv0; n++) {
724     AliESDv0 *v0=GetV0(n);
725     if (v0->GetIndex(0)==last) {
726        v0->SetIndex(0,rm);
727        used--;
728        if (!used) return kTRUE;
729     }
730     if (v0->GetIndex(1)==last) {
731        v0->SetIndex(1,rm);
732        used--;
733        if (!used) return kTRUE;
734     }
735   }
736
737   for (Int_t n=0; n<ncs; n++) {
738     AliESDcascade *cs=GetCascade(n);
739     if (cs->GetIndex()==last) {
740        cs->SetIndex(rm);
741        used--;
742        if (!used) return kTRUE;
743     }
744     AliESDv0 *v0=cs;
745     if (v0->GetIndex(0)==last) {
746        v0->SetIndex(0,rm);
747        used--;
748        if (!used) return kTRUE;
749     }
750     if (v0->GetIndex(1)==last) {
751        v0->SetIndex(1,rm);
752        used--;
753        if (!used) return kTRUE;
754     }
755   }
756
757   for (Int_t n=0; n<nkn; n++) {
758     AliESDkink *kn=GetKink(n);
759     if (kn->GetIndex(0)==last) {
760        kn->SetIndex(rm,0);
761        used--;
762        if (!used) return kTRUE;
763     }
764     if (kn->GetIndex(1)==last) {
765        kn->SetIndex(rm,1);
766        used--;
767        if (!used) return kTRUE;
768     }
769   }
770
771   // Remap the indices of the tracks accosicated with CaloClusters
772   for (Int_t n=0; n<ncl; n++) {
773     AliESDCaloCluster *cluster=GetCaloCluster(n);
774     TArrayI *arr=cluster->GetTracksMatched();
775     Int_t s=arr->GetSize();
776     while (s--) {
777       Int_t idx=arr->At(s);
778       if (idx==last) {
779          arr->AddAt(rm,s);
780          used--; 
781          if (!used) return kTRUE;
782       }
783     }
784   }
785
786   return kTRUE;
787 }
788
789
790 Bool_t AliESDEvent::Clean(Float_t *cleanPars) {
791   //
792   // Remove the data which are not needed for the physics analysis.
793   //
794   // 1) Cleaning the V0 candidates
795   //    ---------------------------
796   //    If the cosine of the V0 pointing angle "csp" and 
797   //    the DCA between the daughter tracks "dca" does not satisfy 
798   //    the conditions 
799   //
800   //     csp > cleanPars[1] + dca/cleanPars[0]*(1.- cleanPars[1])
801   //
802   //    an attempt to remove this V0 candidate from ESD is made.
803   //
804   //    The V0 candidate gets removed if it does not belong to any 
805   //    recosntructed cascade decay
806   //
807   //    12.11.2007, optimal values: cleanPars[0]=0.5, cleanPars[1]=0.999
808   //
809   // 2) Cleaning the tracks
810   //    ----------------------
811   //    If track's transverse parameter is larger than cleanPars[2]
812   //                       OR
813   //    track's longitudinal parameter is larger than cleanPars[3]
814   //    an attempt to remove this track from ESD is made.
815   //
816   //    The track gets removed if it does not come 
817   //    from a reconstructed decay
818   //
819   Bool_t rc=kFALSE;
820
821   Float_t dcaMax=cleanPars[0];
822   Float_t cspMin=cleanPars[1];
823
824   Int_t nV0s=GetNumberOfV0s();
825   for (Int_t i=nV0s-1; i>=0; i--) {
826     AliESDv0 *v0=GetV0(i);
827
828     Float_t dca=v0->GetDcaV0Daughters();
829     Float_t csp=v0->GetV0CosineOfPointingAngle();
830     Float_t cspcut=cspMin + dca/dcaMax*(1.-cspMin);
831     if (csp > cspcut) continue;
832     if (RemoveV0(i)) rc=kTRUE;
833   }
834
835
836   Float_t dmax=cleanPars[2], zmax=cleanPars[3];
837
838   const AliESDVertex *vertex=GetPrimaryVertexSPD();
839   Bool_t vtxOK=vertex->GetStatus();
840   
841   Int_t nTracks=GetNumberOfTracks();
842   for (Int_t i=nTracks-1; i>=0; i--) {
843     AliESDtrack *track=GetTrack(i);
844     Float_t xy,z; track->GetImpactParameters(xy,z);
845     if ((TMath::Abs(xy) > dmax) || (vtxOK && (TMath::Abs(z) > zmax))) {
846       if (RemoveTrack(i)) rc=kTRUE;
847     }
848   }
849
850   return rc;
851 }
852
853 Int_t  AliESDEvent::AddTrack(const AliESDtrack *t) 
854 {
855     // Add track
856     TClonesArray &ftr = *fTracks;
857     AliESDtrack * track = new(ftr[fTracks->GetEntriesFast()])AliESDtrack(*t);
858     track->SetID(fTracks->GetEntriesFast()-1);
859     return  track->GetID();    
860 }
861
862  void AliESDEvent::AddMuonTrack(const AliESDMuonTrack *t) 
863 {
864     TClonesArray &fmu = *fMuonTracks;
865     new(fmu[fMuonTracks->GetEntriesFast()]) AliESDMuonTrack(*t);
866 }
867
868 void AliESDEvent::AddPmdTrack(const AliESDPmdTrack *t) 
869 {
870   TClonesArray &fpmd = *fPmdTracks;
871   new(fpmd[fPmdTracks->GetEntriesFast()]) AliESDPmdTrack(*t);
872 }
873
874 void AliESDEvent::AddTrdTrack(const AliESDTrdTrack *t) 
875 {
876   TClonesArray &ftrd = *fTrdTracks;
877   new(ftrd[fTrdTracks->GetEntriesFast()]) AliESDTrdTrack(*t);
878 }
879
880
881
882
883 Int_t AliESDEvent::AddKink(const AliESDkink *c) 
884 {
885     // Add kink
886     TClonesArray &fk = *fKinks;
887     AliESDkink * kink = new(fk[fKinks->GetEntriesFast()]) AliESDkink(*c);
888     kink->SetID(fKinks->GetEntriesFast()); // CKB different from the other imps..
889     return fKinks->GetEntriesFast()-1;
890 }
891
892
893 void AliESDEvent::AddCascade(const AliESDcascade *c) 
894 {
895   TClonesArray &fc = *fCascades;
896   new(fc[fCascades->GetEntriesFast()]) AliESDcascade(*c);
897 }
898
899
900 Int_t AliESDEvent::AddCaloCluster(const AliESDCaloCluster *c) 
901 {
902     // Add calocluster
903     TClonesArray &fc = *fCaloClusters;
904     AliESDCaloCluster *clus = new(fc[fCaloClusters->GetEntriesFast()]) AliESDCaloCluster(*c);
905     clus->SetID(fCaloClusters->GetEntriesFast()-1);
906     return fCaloClusters->GetEntriesFast()-1;
907   }
908
909
910 void  AliESDEvent::AddRawDataErrorLog(const AliRawDataErrorLog *log) const {
911   TClonesArray &errlogs = *fErrorLogs;
912   new(errlogs[errlogs.GetEntriesFast()])  AliRawDataErrorLog(*log);
913 }
914
915 void  AliESDEvent::SetPrimaryVertexTPC(const AliESDVertex *vertex) 
916 {
917   // Set the TPC vertex
918   // use already allocated space
919   if(fTPCVertex){
920     *fTPCVertex = *vertex;
921     fTPCVertex->SetName(fgkESDListName[kTPCVertex]);
922   }
923 }
924
925 void  AliESDEvent::SetPrimaryVertexSPD(const AliESDVertex *vertex) 
926 {
927   // Set the SPD vertex
928   // use already allocated space
929   if(fSPDVertex){
930     *fSPDVertex = *vertex;
931     fSPDVertex->SetName(fgkESDListName[kSPDVertex]);
932   }
933 }
934
935 void  AliESDEvent::SetPrimaryVertexTracks(const AliESDVertex *vertex) 
936 {
937   // Set the primary vertex reconstructed using he ESD tracks.
938   // use already allocated space
939   if(fPrimaryVertex){
940     *fPrimaryVertex = *vertex;
941     fPrimaryVertex->SetName(fgkESDListName[kPrimaryVertex]);
942   }
943 }
944
945 const AliESDVertex * AliESDEvent::GetPrimaryVertex() const 
946 {
947   //
948   // Get the "best" available reconstructed primary vertex.
949   //
950   if(fPrimaryVertex){
951     if (fPrimaryVertex->GetStatus()) return fPrimaryVertex;
952   }
953   if(fSPDVertex){
954     if (fSPDVertex->GetStatus()) return fSPDVertex;
955   }
956   if(fTPCVertex) return fTPCVertex;
957   
958   AliWarning("No primary vertex available. Returning the \"default\"...");
959   return fSPDVertex;
960 }
961
962 void AliESDEvent::SetMultiplicity(const AliMultiplicity *mul) 
963 {
964   // Set the SPD Multiplicity
965   if(fSPDMult){
966     *fSPDMult = *mul;
967   }
968 }
969
970
971 void AliESDEvent::SetFMDData(AliESDFMD * obj) 
972
973   // use already allocated space
974   if(fESDFMD){
975     *fESDFMD = *obj;
976   }
977 }
978
979 void AliESDEvent::SetVZEROData(AliESDVZERO * obj)
980
981   // use already allocated space
982   if(fESDVZERO)
983     *fESDVZERO = *obj;
984 }
985
986 void AliESDEvent::SetACORDEData(AliESDACORDE * obj)
987 {
988   if(fESDACORDE)
989     *fESDACORDE = *obj;
990 }
991
992
993 void AliESDEvent::GetESDfriend(AliESDfriend *ev) const 
994 {
995   //
996   // Extracts the complementary info from the ESD
997   //
998   if (!ev) return;
999
1000   Int_t ntrk=GetNumberOfTracks();
1001
1002   for (Int_t i=0; i<ntrk; i++) {
1003     AliESDtrack *t=GetTrack(i);
1004     const AliESDfriendTrack *f=t->GetFriendTrack();
1005     ev->AddTrack(f);
1006
1007     t->ReleaseESDfriendTrack();// Not to have two copies of "friendTrack"
1008
1009   }
1010
1011   AliESDfriend *fr = (AliESDfriend*)(const_cast<AliESDEvent*>(this)->FindListObject("AliESDfriend"));
1012   if (fr) ev->SetVZEROfriend(fr->GetVZEROfriend());
1013 }
1014
1015 void AliESDEvent::AddObject(TObject* obj) 
1016 {
1017   // Add an object to the list of object.
1018   // Please be aware that in order to increase performance you should
1019   // refrain from using TObjArrays (if possible). Use TClonesArrays, instead.
1020   fESDObjects->SetOwner(kTRUE);
1021   fESDObjects->AddLast(obj);
1022 }
1023
1024
1025 void AliESDEvent::GetStdContent() 
1026 {
1027   // set pointers for standard content
1028   // get by name much safer and not a big overhead since not called very often
1029  
1030   fESDRun = (AliESDRun*)fESDObjects->FindObject(fgkESDListName[kESDRun]);
1031   fHeader = (AliESDHeader*)fESDObjects->FindObject(fgkESDListName[kHeader]);
1032   fESDZDC = (AliESDZDC*)fESDObjects->FindObject(fgkESDListName[kESDZDC]);
1033   fESDFMD = (AliESDFMD*)fESDObjects->FindObject(fgkESDListName[kESDFMD]);
1034   fESDVZERO = (AliESDVZERO*)fESDObjects->FindObject(fgkESDListName[kESDVZERO]);
1035   fESDTZERO = (AliESDTZERO*)fESDObjects->FindObject(fgkESDListName[kESDTZERO]);
1036   fTPCVertex = (AliESDVertex*)fESDObjects->FindObject(fgkESDListName[kTPCVertex]);
1037   fSPDVertex = (AliESDVertex*)fESDObjects->FindObject(fgkESDListName[kSPDVertex]);
1038   fPrimaryVertex = (AliESDVertex*)fESDObjects->FindObject(fgkESDListName[kPrimaryVertex]);
1039   fSPDMult =       (AliMultiplicity*)fESDObjects->FindObject(fgkESDListName[kSPDMult]);
1040   fPHOSTrigger = (AliESDCaloTrigger*)fESDObjects->FindObject(fgkESDListName[kPHOSTrigger]);
1041   fEMCALTrigger = (AliESDCaloTrigger*)fESDObjects->FindObject(fgkESDListName[kEMCALTrigger]);
1042   fTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kTracks]);
1043   fMuonTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kMuonTracks]);
1044   fPmdTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kPmdTracks]);
1045   fTrdTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kTrdTracks]);
1046   fV0s = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kV0s]);
1047   fCascades = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kCascades]);
1048   fKinks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kKinks]);
1049   fCaloClusters = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kCaloClusters]);
1050   fEMCALCells = (AliESDCaloCells*)fESDObjects->FindObject(fgkESDListName[kEMCALCells]);
1051   fPHOSCells = (AliESDCaloCells*)fESDObjects->FindObject(fgkESDListName[kPHOSCells]);
1052   fErrorLogs = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kErrorLogs]);
1053   fESDACORDE = (AliESDACORDE*)fESDObjects->FindObject(fgkESDListName[kESDACORDE]);
1054
1055 }
1056
1057 void AliESDEvent::SetStdNames(){
1058   // Set the names of the standard contents
1059   // 
1060   if(fESDObjects->GetEntries()>=kESDListN){
1061     for(int i = 0;i < fESDObjects->GetEntries() && i<kESDListN;i++){
1062       TObject *fObj = fESDObjects->At(i);
1063       if(fObj->InheritsFrom("TNamed")){
1064         ((TNamed*)fObj)->SetName(fgkESDListName[i]);
1065       }
1066       else if(fObj->InheritsFrom("TClonesArray")){
1067         ((TClonesArray*)fObj)->SetName(fgkESDListName[i]);
1068       }
1069     }
1070   }
1071   else{
1072      AliWarning("Std Entries missing");
1073   }
1074
1075
1076
1077 void AliESDEvent::CreateStdContent(Bool_t bUseThisList){
1078   fUseOwnList = bUseThisList;
1079   CreateStdContent();
1080 }
1081
1082 void AliESDEvent::CreateStdContent() 
1083 {
1084   // create the standard AOD content and set pointers
1085
1086   // create standard objects and add them to the TList of objects
1087   AddObject(new AliESDRun());
1088   AddObject(new AliESDHeader());
1089   AddObject(new AliESDZDC());
1090   AddObject(new AliESDFMD());
1091   AddObject(new AliESDVZERO());
1092   AddObject(new AliESDTZERO());
1093   AddObject(new AliESDVertex());
1094   AddObject(new AliESDVertex());
1095   AddObject(new AliESDVertex());
1096   AddObject(new AliMultiplicity());
1097   AddObject(new AliESDCaloTrigger());
1098   AddObject(new AliESDCaloTrigger());
1099   AddObject(new TClonesArray("AliESDtrack",0));
1100   AddObject(new TClonesArray("AliESDMuonTrack",0));
1101   AddObject(new TClonesArray("AliESDPmdTrack",0));
1102   AddObject(new TClonesArray("AliESDTrdTrack",0));
1103   AddObject(new TClonesArray("AliESDv0",0));
1104   AddObject(new TClonesArray("AliESDcascade",0));
1105   AddObject(new TClonesArray("AliESDkink",0));
1106   AddObject(new TClonesArray("AliESDCaloCluster",0));
1107   AddObject(new AliESDCaloCells());
1108   AddObject(new AliESDCaloCells());
1109   AddObject(new TClonesArray("AliRawDataErrorLog",0));
1110   AddObject(new AliESDACORDE()); 
1111
1112   // check the order of the indices against enum...
1113
1114   // set names
1115   SetStdNames();
1116   // read back pointers
1117   GetStdContent();
1118 }
1119
1120 TObject* AliESDEvent::FindListObject(const char *name){
1121 //
1122 // Find object with name "name" in the list of branches
1123 //
1124   if(fESDObjects){
1125     return fESDObjects->FindObject(name);
1126   }
1127   return 0;
1128
1129
1130 Int_t AliESDEvent::GetPHOSClusters(TRefArray *clusters) const
1131 {
1132   // fills the provided TRefArray with all found phos clusters
1133   
1134   clusters->Clear();
1135   
1136   AliESDCaloCluster *cl = 0;
1137   for (Int_t i = 0; i < GetNumberOfCaloClusters(); i++) {
1138     
1139     if ( (cl = GetCaloCluster(i)) ) {
1140       if (cl->IsPHOS()){
1141         clusters->Add(cl);
1142         AliDebug(1,Form("IsPHOS cluster %d Size: %d \n",i,clusters->GetEntriesFast()));
1143       }
1144     }
1145   }
1146   return clusters->GetEntriesFast();
1147 }
1148
1149 Int_t AliESDEvent::GetEMCALClusters(TRefArray *clusters) const
1150 {
1151   // fills the provided TRefArray with all found emcal clusters
1152
1153   clusters->Clear();
1154
1155   AliESDCaloCluster *cl = 0;
1156   for (Int_t i = 0; i < GetNumberOfCaloClusters(); i++) {
1157
1158     if ( (cl = GetCaloCluster(i)) ) {
1159       if (cl->IsEMCAL()){
1160         clusters->Add(cl);
1161         AliDebug(1,Form("IsEMCAL cluster %d Size: %d \n",i,clusters->GetEntriesFast()));
1162       }
1163     }
1164   }
1165   return clusters->GetEntriesFast();
1166 }
1167
1168 void AliESDEvent::WriteToTree(TTree* tree) const {
1169   // Book the branches as in TTree::Branch(TCollection*)
1170   // but add a "." at the end of top level branches which are
1171   // not a TClonesArray
1172
1173
1174   TString branchname;
1175   TIter next(fESDObjects);
1176   const Int_t kSplitlevel = 99; // default value in TTree::Branch()
1177   const Int_t kBufsize = 32000; // default value in TTree::Branch()
1178   TObject *obj = 0;
1179
1180   while ((obj = next())) {
1181     branchname.Form("%s", obj->GetName());
1182     if(branchname.CompareTo("AliESDfriend")==0)branchname = "ESDfriend.";
1183     if ((kSplitlevel > 1) &&  !obj->InheritsFrom(TClonesArray::Class())) {
1184       if(!branchname.EndsWith("."))branchname += ".";
1185     }
1186     if (!tree->FindBranch(branchname)) {
1187       tree->Bronch(branchname, obj->ClassName(), fESDObjects->GetObjectRef(obj),
1188                    kBufsize, kSplitlevel - 1);
1189     }
1190   }
1191 }
1192
1193
1194 void AliESDEvent::ReadFromTree(TTree *tree, Option_t* /*opt*/){
1195 //
1196 // Connect the ESDEvent to a tree
1197 //
1198   if(!tree){
1199     AliWarning("AliESDEvent::ReadFromTree() Zero Pointer to Tree \n");
1200     return;
1201   }
1202   // load the TTree
1203   if(!tree->GetTree())tree->LoadTree(0);
1204
1205   // if we find the "ESD" branch on the tree we do have the old structure
1206   if(tree->GetBranch("ESD")) {
1207     char ** address  = (char **)(tree->GetBranch("ESD")->GetAddress());
1208     // do we have the friend branch
1209     TBranch * esdFB = tree->GetBranch("ESDfriend.");
1210     char ** addressF = 0;
1211     if(esdFB)addressF = (char **)(esdFB->GetAddress());
1212     if (!address) {
1213       AliInfo("AliESDEvent::ReadFromTree() Reading old Tree");
1214       tree->SetBranchAddress("ESD",       &fESDOld);
1215       if(esdFB){
1216         tree->SetBranchAddress("ESDfriend.",&fESDFriendOld);
1217       }
1218     } else {
1219       AliInfo("AliESDEvent::ReadFromTree() Reading old Tree");
1220       AliInfo("Branch already connected. Using existing branch address.");
1221       fESDOld       = (AliESD*)       (*address);
1222       // addressF can still be 0, since branch needs to switched on
1223       if(addressF)fESDFriendOld = (AliESDfriend*) (*addressF);
1224     }
1225                                        
1226     //  have already connected the old ESD structure... ?
1227     // reuse also the pointer of the AlliESDEvent
1228     // otherwise create new ones
1229     TList* connectedList = (TList*) (tree->GetUserInfo()->FindObject("ESDObjectsConnectedToTree"));
1230   
1231     if(connectedList){
1232       // If connected use the connected list of objects
1233       if(fESDObjects!= connectedList){
1234         // protect when called twice 
1235         fESDObjects->Delete();
1236         fESDObjects = connectedList;
1237       }
1238       GetStdContent(); 
1239
1240       
1241       // The pointer to the friend changes when called twice via InitIO
1242       // since AliESDEvent is deleted
1243       TObject* oldf = FindListObject("AliESDfriend");
1244       TObject* newf = 0;
1245       if(addressF){
1246         newf = (TObject*)*addressF;
1247       }
1248       if(newf!=0&&oldf!=newf){
1249         // remove the old reference
1250         // Should we also delete it? Or is this handled in TTree I/O
1251         // since it is created by the first SetBranchAddress
1252         fESDObjects->Remove(oldf);
1253         // add the new one 
1254         fESDObjects->Add(newf);
1255       }
1256       
1257       fConnected = true;
1258       return;
1259     }
1260     // else...    
1261     CreateStdContent(); // create for copy
1262     // if we have the esdfriend add it, so we always can access it via the userinfo
1263     if(fESDFriendOld)AddObject(fESDFriendOld);
1264     // we are not owner of the list objects 
1265     // must not delete it
1266     fESDObjects->SetOwner(kFALSE);
1267     fESDObjects->SetName("ESDObjectsConnectedToTree");
1268     tree->GetUserInfo()->Add(fESDObjects);
1269     fConnected = true;
1270     return;
1271   }
1272   
1273
1274     delete fESDOld;
1275     fESDOld = 0;
1276   // Try to find AliESDEvent
1277   AliESDEvent *esdEvent = 0;
1278   esdEvent = (AliESDEvent*)tree->GetTree()->GetUserInfo()->FindObject("AliESDEvent");
1279   if(esdEvent){   
1280       // Check if already connected to tree
1281     TList* connectedList = (TList*) (tree->GetUserInfo()->FindObject("ESDObjectsConnectedToTree"));
1282     if (connectedList) {
1283       // If connected use the connected list if objects
1284       fESDObjects->Delete();
1285       fESDObjects = connectedList;
1286       GetStdContent(); 
1287       fConnected = true;
1288       return;
1289     }
1290
1291     // Connect to tree
1292     // prevent a memory leak when reading back the TList
1293
1294     if(!fUseOwnList){
1295       delete fESDObjects;
1296       fESDObjects = 0;
1297       // create a new TList from the UserInfo TList... 
1298       // copy constructor does not work...
1299       fESDObjects = (TList*)(esdEvent->GetList()->Clone());
1300       fESDObjects->SetOwner(kFALSE);
1301     }
1302     else if ( fESDObjects->GetEntries()==0){
1303       // at least create the std content if we want to read to our list
1304       CreateStdContent(); 
1305     }
1306
1307     // in principle
1308     // we only need new things in the list if we do no already have it..
1309     // TODO just add new entries
1310
1311     if(fESDObjects->GetEntries()<kESDListN){
1312       AliWarning(Form("AliESDEvent::ReadFromTree() TList contains less than the standard contents %d < %d \n",
1313                       fESDObjects->GetEntries(),kESDListN));
1314     }
1315     // set the branch addresses
1316     TIter next(fESDObjects);
1317     TNamed *el;
1318     while((el=(TNamed*)next())){
1319       TString bname(el->GetName());
1320       if(bname.CompareTo("AliESDfriend")==0)
1321         {
1322           // AliESDfriend does not have a name ...
1323           tree->SetBranchAddress("ESDfriend.",fESDObjects->GetObjectRef(el));
1324         }
1325       else{
1326         // check if branch exists under this Name
1327         TBranch *br = tree->GetBranch(bname.Data());
1328         if(br){
1329           tree->SetBranchAddress(bname.Data(),fESDObjects->GetObjectRef(el));
1330         }
1331         else{
1332           br = tree->GetBranch(Form("%s.",bname.Data()));
1333           if(br){
1334             tree->SetBranchAddress(Form("%s.",bname.Data()),fESDObjects->GetObjectRef(el));
1335           }
1336           else{
1337             AliWarning(Form("AliESDEvent::ReadFromTree() No Branch found with Name %s or %s.",bname.Data(),bname.Data()));
1338           }
1339
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     fESDObjects->SetName("ESDObjectsConnectedToTree");
1348     // we are not owner of the list objects 
1349     // must not delete it
1350     tree->GetUserInfo()->Add(fESDObjects);
1351     fConnected = true;
1352   }// no esdEvent
1353   else {
1354     // we can't get the list from the user data, create standard content
1355     // and set it by hand (no ESDfriend at the moment
1356     CreateStdContent();
1357     TIter next(fESDObjects);
1358     TNamed *el;
1359     while((el=(TNamed*)next())){
1360       TString bname(el->GetName());    
1361       TBranch *br = tree->GetBranch(bname.Data());
1362       if(br){
1363         tree->SetBranchAddress(bname.Data(),fESDObjects->GetObjectRef(el));
1364       }
1365       else{
1366         br = tree->GetBranch(Form("%s.",bname.Data()));
1367         if(br){
1368           tree->SetBranchAddress(Form("%s.",bname.Data()),fESDObjects->GetObjectRef(el));
1369         }
1370       }
1371     }
1372     GetStdContent();
1373     // when reading back we are not owner of the list 
1374     // must not delete it
1375     fESDObjects->SetOwner(kFALSE);
1376   }
1377 }
1378
1379
1380 void AliESDEvent::CopyFromOldESD()
1381 {
1382   // Method which copies over everthing from the old esd structure to the 
1383   // new  
1384   if(fESDOld){
1385     ResetStdContent();
1386      // Run
1387     SetRunNumber(fESDOld->GetRunNumber());
1388     SetPeriodNumber(fESDOld->GetPeriodNumber());
1389     SetMagneticField(fESDOld->GetMagneticField());
1390   
1391     // leave out diamond ...
1392     // SetDiamond(const AliESDVertex *vertex) { fESDRun->SetDiamond(vertex);}
1393
1394     // header
1395     SetTriggerMask(fESDOld->GetTriggerMask());
1396     SetOrbitNumber(fESDOld->GetOrbitNumber());
1397     SetTimeStamp(fESDOld->GetTimeStamp());
1398     SetEventType(fESDOld->GetEventType());
1399     SetEventNumberInFile(fESDOld->GetEventNumberInFile());
1400     SetBunchCrossNumber(fESDOld->GetBunchCrossNumber());
1401     SetTriggerCluster(fESDOld->GetTriggerCluster());
1402
1403     // ZDC
1404
1405     SetZDC(fESDOld->GetZDCN1Energy(),
1406            fESDOld->GetZDCP1Energy(),
1407            fESDOld->GetZDCEMEnergy(),
1408            0,
1409            fESDOld->GetZDCN2Energy(),
1410            fESDOld->GetZDCP2Energy(),
1411            fESDOld->GetZDCParticipants(),
1412            0);
1413
1414     // FMD
1415     
1416     if(fESDOld->GetFMDData())SetFMDData(fESDOld->GetFMDData());
1417
1418     // T0
1419
1420     SetT0zVertex(fESDOld->GetT0zVertex());
1421     SetT0(fESDOld->GetT0());
1422     //  leave amps out
1423
1424     // VZERO
1425     if (fESDOld->GetVZEROData()) SetVZEROData(fESDOld->GetVZEROData());
1426
1427     if(fESDOld->GetVertex())SetPrimaryVertexSPD(fESDOld->GetVertex());
1428
1429     if(fESDOld->GetPrimaryVertex())SetPrimaryVertexTracks(fESDOld->GetPrimaryVertex());
1430
1431     if(fESDOld->GetMultiplicity())SetMultiplicity(fESDOld->GetMultiplicity());
1432
1433     for(int i = 0;i<fESDOld->GetNumberOfTracks();i++){
1434       AddTrack(fESDOld->GetTrack(i));
1435     }
1436
1437     for(int i = 0;i<fESDOld->GetNumberOfMuonTracks();i++){
1438       AddMuonTrack(fESDOld->GetMuonTrack(i));
1439     }
1440
1441     for(int i = 0;i<fESDOld->GetNumberOfPmdTracks();i++){
1442       AddPmdTrack(fESDOld->GetPmdTrack(i));
1443     }
1444
1445     for(int i = 0;i<fESDOld->GetNumberOfTrdTracks();i++){
1446       AddTrdTrack(fESDOld->GetTrdTrack(i));
1447     }
1448
1449     for(int i = 0;i<fESDOld->GetNumberOfV0s();i++){
1450       AddV0(fESDOld->GetV0(i));
1451     }
1452
1453     for(int i = 0;i<fESDOld->GetNumberOfCascades();i++){
1454       AddCascade(fESDOld->GetCascade(i));
1455     }
1456
1457     for(int i = 0;i<fESDOld->GetNumberOfKinks();i++){
1458       AddKink(fESDOld->GetKink(i));
1459     }
1460
1461
1462     for(int i = 0;i<fESDOld->GetNumberOfCaloClusters();i++){
1463       AddCaloCluster(fESDOld->GetCaloCluster(i));
1464     }
1465
1466   }// if fesdold
1467 }
1468
1469
1470