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