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