]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliESDEvent.cxx
f8d8026f8dcb90165c3d06b1beba2082fbd0fead
[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 a reconstructed decay
530   Int_t nv0=GetNumberOfV0s();
531   for (Int_t n=0; n<nv0; n++) {
532     AliESDv0 *v0=GetV0(n);
533
534     Int_t idx=v0->GetNindex();
535     if (rm==idx) return kFALSE;
536     if (idx==last) used++;
537
538     idx=v0->GetPindex();
539     if (rm==idx) return kFALSE;
540     if (idx==last) used++;
541   }
542
543   Int_t ncs=GetNumberOfCascades();
544   for (Int_t n=0; n<ncs; n++) {
545     AliESDcascade *cs=GetCascade(n);
546
547     Int_t idx=cs->GetIndex();
548     if (rm==idx) return kFALSE;
549     if (idx==last) used++;
550   }
551
552   Int_t nkn=GetNumberOfKinks();
553   for (Int_t n=0; n<nkn; n++) {
554     AliESDkink *kn=GetKink(n);
555
556     Int_t idx=kn->GetIndex(0);
557     if (rm==idx) return kFALSE;
558     if (idx==last) used++;
559
560     idx=kn->GetIndex(1);
561     if (rm==idx) return kFALSE;
562     if (idx==last) used++;
563   }
564
565
566   //Replace the removed track with the last track 
567   TClonesArray &a=*fTracks;
568   delete a.RemoveAt(rm);
569
570   if (rm==last) return kTRUE;
571
572   AliESDtrack *t=GetTrack(last);
573   t->SetID(rm);
574   new (a[rm]) AliESDtrack(*t);
575   delete a.RemoveAt(last);
576
577   if (!used) return kTRUE;
578   
579
580   // Remap the indices of the daughters of reconstructed decays
581   for (Int_t n=0; n<nv0; n++) {
582     AliESDv0 *v0=GetV0(n);
583     if (v0->GetIndex(0)==last) {
584        v0->SetIndex(0,rm);
585        used--;
586        if (!used) return kTRUE;
587     }
588     if (v0->GetIndex(1)==last) {
589        v0->SetIndex(1,rm);
590        used--;
591        if (!used) return kTRUE;
592     }
593   }
594
595   for (Int_t n=0; n<ncs; n++) {
596     AliESDcascade *cs=GetCascade(n);
597     if (cs->GetIndex()==last) {
598        cs->SetIndex(rm);
599        used--;
600        if (!used) return kTRUE;
601     }
602   }
603
604   for (Int_t n=0; n<nkn; n++) {
605     AliESDkink *kn=GetKink(n);
606     if (kn->GetIndex(0)==last) {
607        kn->SetIndex(rm,0);
608        used--;
609        if (!used) return kTRUE;
610     }
611     if (kn->GetIndex(1)==last) {
612        kn->SetIndex(rm,1);
613        used--;
614        if (!used) return kTRUE;
615     }
616   }
617
618   return kTRUE;
619 }
620
621
622 Bool_t AliESDEvent::Clean(Float_t *cleanPars) {
623   //
624   // Remove the data which are not needed for the physics analysis.
625   //
626   // 1) Cleaning the V0 candidates
627   //    ---------------------------
628   //    If the cosine of the V0 pointing angle "csp" and 
629   //    the DCA between the daughter tracks "dca" does not satisfy 
630   //    the conditions 
631   //
632   //     csp > cleanPars[1] + dca/cleanPars[0]*(1.- cleanPars[1])
633   //
634   //    an attempt to remove this V0 candidate from ESD is made.
635   //
636   //    The V0 candidate gets removed if it does not belong to any 
637   //    recosntructed cascade decay
638   //
639   //    12.11.2007, optimal values: cleanPars[0]=0.5, cleanPars[1]=0.999
640   //
641   // 2) Cleaning the tracks
642   //    ----------------------
643   //    If track's transverse parameter is larger than cleanPars[2]
644   //                       OR
645   //    track's longitudinal parameter is larger than cleanPars[3]
646   //    an attempt to remove this track from ESD is made.
647   //
648   //    The track gets removed if it does not come 
649   //    from a reconstructed decay
650   //
651   Bool_t rc=kFALSE;
652
653   Float_t dcaMax=cleanPars[0];
654   Float_t cspMin=cleanPars[1];
655
656   Int_t nV0s=GetNumberOfV0s();
657   for (Int_t i=nV0s-1; i>=0; i--) {
658     AliESDv0 *v0=GetV0(i);
659
660     Float_t dca=v0->GetDcaV0Daughters();
661     Float_t csp=v0->GetV0CosineOfPointingAngle();
662     Float_t cspcut=cspMin + dca/dcaMax*(1.-cspMin);
663     if (csp > cspcut) continue;
664     if (RemoveV0(i)) rc=kTRUE;
665   }
666
667
668   Float_t dmax=cleanPars[2], zmax=cleanPars[3];
669
670   const AliESDVertex *vertex=GetPrimaryVertexSPD();
671   Bool_t vtxOK=vertex->GetStatus();
672   
673   Int_t nTracks=GetNumberOfTracks();
674   for (Int_t i=nTracks-1; i>=0; i--) {
675     AliESDtrack *track=GetTrack(i);
676     Float_t xy,z; track->GetImpactParameters(xy,z);
677     if ((TMath::Abs(xy) > dmax) || (vtxOK && (TMath::Abs(z) > zmax))) {
678       if (RemoveTrack(i)) rc=kTRUE;
679     }
680   }
681
682   return rc;
683 }
684
685 Int_t  AliESDEvent::AddTrack(const AliESDtrack *t) 
686 {
687     // Add track
688     TClonesArray &ftr = *fTracks;
689     AliESDtrack * track = new(ftr[fTracks->GetEntriesFast()])AliESDtrack(*t);
690     track->SetID(fTracks->GetEntriesFast()-1);
691     return  track->GetID();    
692 }
693
694  void AliESDEvent::AddMuonTrack(const AliESDMuonTrack *t) 
695 {
696     TClonesArray &fmu = *fMuonTracks;
697     new(fmu[fMuonTracks->GetEntriesFast()]) AliESDMuonTrack(*t);
698 }
699
700 void AliESDEvent::AddPmdTrack(const AliESDPmdTrack *t) 
701 {
702   TClonesArray &fpmd = *fPmdTracks;
703   new(fpmd[fPmdTracks->GetEntriesFast()]) AliESDPmdTrack(*t);
704 }
705
706 void AliESDEvent::AddTrdTrack(const AliESDTrdTrack *t) 
707 {
708   TClonesArray &ftrd = *fTrdTracks;
709   new(ftrd[fTrdTracks->GetEntriesFast()]) AliESDTrdTrack(*t);
710 }
711
712
713
714
715 Int_t AliESDEvent::AddKink(const AliESDkink *c) 
716 {
717     // Add kink
718     TClonesArray &fk = *fKinks;
719     AliESDkink * kink = new(fk[fKinks->GetEntriesFast()]) AliESDkink(*c);
720     kink->SetID(fKinks->GetEntriesFast()); // CKB different from the other imps..
721     return fKinks->GetEntriesFast()-1;
722 }
723
724
725 void AliESDEvent::AddCascade(const AliESDcascade *c) 
726 {
727   TClonesArray &fc = *fCascades;
728   new(fc[fCascades->GetEntriesFast()]) AliESDcascade(*c);
729 }
730
731
732 Int_t AliESDEvent::AddCaloCluster(const AliESDCaloCluster *c) 
733 {
734     // Add calocluster
735     TClonesArray &fc = *fCaloClusters;
736     AliESDCaloCluster *clus = new(fc[fCaloClusters->GetEntriesFast()]) AliESDCaloCluster(*c);
737     clus->SetID(fCaloClusters->GetEntriesFast()-1);
738     return fCaloClusters->GetEntriesFast()-1;
739   }
740
741
742 void  AliESDEvent::AddRawDataErrorLog(const AliRawDataErrorLog *log) const {
743   TClonesArray &errlogs = *fErrorLogs;
744   new(errlogs[errlogs.GetEntriesFast()])  AliRawDataErrorLog(*log);
745 }
746
747 void  AliESDEvent::SetPrimaryVertexTPC(const AliESDVertex *vertex) 
748 {
749   // Set the TPC vertex
750   // use already allocated space
751   if(fTPCVertex){
752     *fTPCVertex = *vertex;
753     fTPCVertex->SetName(fgkESDListName[kTPCVertex]);
754   }
755 }
756
757 void  AliESDEvent::SetPrimaryVertexSPD(const AliESDVertex *vertex) 
758 {
759   // Set the SPD vertex
760   // use already allocated space
761   if(fSPDVertex){
762     *fSPDVertex = *vertex;
763     fSPDVertex->SetName(fgkESDListName[kSPDVertex]);
764   }
765 }
766
767 void  AliESDEvent::SetPrimaryVertex(const AliESDVertex *vertex) 
768 {
769   // Set the primary vertex
770   // use already allocated space
771   if(fPrimaryVertex){
772     *fPrimaryVertex = *vertex;
773     fPrimaryVertex->SetName(fgkESDListName[kPrimaryVertex]);
774   }
775 }
776
777 void AliESDEvent::SetMultiplicity(const AliMultiplicity *mul) 
778 {
779   // Set the SPD Multiplicity
780   if(fSPDMult){
781     *fSPDMult = *mul;
782   }
783 }
784
785
786 void AliESDEvent::SetFMDData(AliESDFMD * obj) 
787
788   // use already allocated space
789   if(fESDFMD){
790     *fESDFMD = *obj;
791   }
792 }
793
794 void AliESDEvent::SetVZEROData(AliESDVZERO * obj)
795
796   // use already allocated space
797   if(fESDVZERO)
798     *fESDVZERO = *obj;
799 }
800
801 void AliESDEvent::GetESDfriend(AliESDfriend *ev) const 
802 {
803   //
804   // Extracts the complementary info from the ESD
805   //
806   if (!ev) return;
807
808   Int_t ntrk=GetNumberOfTracks();
809
810   for (Int_t i=0; i<ntrk; i++) {
811     AliESDtrack *t=GetTrack(i);
812     const AliESDfriendTrack *f=t->GetFriendTrack();
813     ev->AddTrack(f);
814
815     t->ReleaseESDfriendTrack();// Not to have two copies of "friendTrack"
816
817   }
818 }
819
820
821 void AliESDEvent::AddObject(TObject* obj) 
822 {
823   // Add an object to the list of object.
824   // Please be aware that in order to increase performance you should
825   // refrain from using TObjArrays (if possible). Use TClonesArrays, instead.
826   fESDObjects->SetOwner(kTRUE);
827   fESDObjects->AddLast(obj);
828 }
829
830
831 void AliESDEvent::GetStdContent() 
832 {
833   // set pointers for standard content
834   // get by name much safer and not a big overhead since not called very often
835  
836   fESDRun = (AliESDRun*)fESDObjects->FindObject(fgkESDListName[kESDRun]);
837   fHeader = (AliESDHeader*)fESDObjects->FindObject(fgkESDListName[kHeader]);
838   fESDZDC = (AliESDZDC*)fESDObjects->FindObject(fgkESDListName[kESDZDC]);
839   fESDFMD = (AliESDFMD*)fESDObjects->FindObject(fgkESDListName[kESDFMD]);
840   fESDVZERO = (AliESDVZERO*)fESDObjects->FindObject(fgkESDListName[kESDVZERO]);
841   fESDTZERO = (AliESDTZERO*)fESDObjects->FindObject(fgkESDListName[kESDTZERO]);
842   fTPCVertex = (AliESDVertex*)fESDObjects->FindObject(fgkESDListName[kTPCVertex]);
843   fSPDVertex = (AliESDVertex*)fESDObjects->FindObject(fgkESDListName[kSPDVertex]);
844   fPrimaryVertex = (AliESDVertex*)fESDObjects->FindObject(fgkESDListName[kPrimaryVertex]);
845   fSPDMult =       (AliMultiplicity*)fESDObjects->FindObject(fgkESDListName[kSPDMult]);
846   fPHOSTrigger = (AliESDCaloTrigger*)fESDObjects->FindObject(fgkESDListName[kPHOSTrigger]);
847   fEMCALTrigger = (AliESDCaloTrigger*)fESDObjects->FindObject(fgkESDListName[kEMCALTrigger]);
848   fTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kTracks]);
849   fMuonTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kMuonTracks]);
850   fPmdTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kPmdTracks]);
851   fTrdTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kTrdTracks]);
852   fV0s = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kV0s]);
853   fCascades = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kCascades]);
854   fKinks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kKinks]);
855   fCaloClusters = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kCaloClusters]);
856   fEMCALCells = (AliESDCaloCells*)fESDObjects->FindObject(fgkESDListName[kEMCALCells]);
857   fPHOSCells = (AliESDCaloCells*)fESDObjects->FindObject(fgkESDListName[kPHOSCells]);
858   fErrorLogs = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kErrorLogs]);
859
860 }
861
862 void AliESDEvent::SetStdNames(){
863   // Set the names of the standard contents
864   // 
865   if(fESDObjects->GetEntries()==kESDListN){
866     for(int i = 0;i < fESDObjects->GetEntries();i++){
867       TObject *fObj = fESDObjects->At(i);
868       if(fObj->InheritsFrom("TNamed")){
869         ((TNamed*)fObj)->SetName(fgkESDListName[i]);
870       }
871       else if(fObj->InheritsFrom("TClonesArray")){
872         ((TClonesArray*)fObj)->SetName(fgkESDListName[i]);
873       }
874     }
875   }
876   else{
877     printf("%s:%d SetStdNames() Wrong number of Std Entries \n",(char*)__FILE__,__LINE__);
878   }
879
880
881 void AliESDEvent::CreateStdContent() 
882 {
883   // create the standard AOD content and set pointers
884
885   // create standard objects and add them to the TList of objects
886   AddObject(new AliESDRun());
887   AddObject(new AliESDHeader());
888   AddObject(new AliESDZDC());
889   AddObject(new AliESDFMD());
890   AddObject(new AliESDVZERO());
891   AddObject(new AliESDTZERO());
892   AddObject(new AliESDVertex());
893   AddObject(new AliESDVertex());
894   AddObject(new AliESDVertex());
895   AddObject(new AliMultiplicity());
896   AddObject(new AliESDCaloTrigger());
897   AddObject(new AliESDCaloTrigger());
898   AddObject(new TClonesArray("AliESDtrack",0));
899   AddObject(new TClonesArray("AliESDMuonTrack",0));
900   AddObject(new TClonesArray("AliESDPmdTrack",0));
901   AddObject(new TClonesArray("AliESDTrdTrack",0));
902   AddObject(new TClonesArray("AliESDv0",0));
903   AddObject(new TClonesArray("AliESDcascade",0));
904   AddObject(new TClonesArray("AliESDkink",0));
905   AddObject(new TClonesArray("AliESDCaloCluster",0));
906   AddObject(new AliESDCaloCells());
907   AddObject(new AliESDCaloCells());
908   AddObject(new TClonesArray("AliRawDataErrorLog",0));
909
910   // check the order of the indices against enum...
911
912   // set names
913   SetStdNames();
914   // read back pointers
915   GetStdContent();
916 }
917
918 TObject* AliESDEvent::FindListObject(const char *name){
919 //
920 // Find object with name "name" in the list of branches
921 //
922   if(fESDObjects){
923     return fESDObjects->FindObject(name);
924   }
925   return 0;
926
927
928 Int_t AliESDEvent::GetPHOSClusters(TRefArray *clusters) const
929 {
930   // fills the provided TRefArray with all found phos clusters
931   
932   clusters->Clear();
933   
934   AliESDCaloCluster *cl = 0;
935   for (Int_t i = 0; i < GetNumberOfCaloClusters(); i++) {
936     
937     if ( (cl = GetCaloCluster(i)) ) {
938       if (cl->IsPHOS()){
939         clusters->Add(cl);
940         AliDebug(1,Form("IsPHOS cluster %d Size: %d \n",i,clusters->GetEntriesFast()));
941       }
942     }
943   }
944   return clusters->GetEntriesFast();
945 }
946
947 Int_t AliESDEvent::GetEMCALClusters(TRefArray *clusters) const
948 {
949   // fills the provided TRefArray with all found emcal clusters
950
951   clusters->Clear();
952
953   AliESDCaloCluster *cl = 0;
954   for (Int_t i = 0; i < GetNumberOfCaloClusters(); i++) {
955
956     if ( (cl = GetCaloCluster(i)) ) {
957       if (cl->IsEMCAL()){
958         clusters->Add(cl);
959         AliDebug(1,Form("IsEMCAL cluster %d Size: %d \n",i,clusters->GetEntriesFast()));
960       }
961     }
962   }
963   return clusters->GetEntriesFast();
964 }
965
966 const void AliESDEvent::WriteToTree(TTree* tree) const {
967   // Book the branches as in TTree::Branch(TCollection*)
968   // but add a "." at the end of top level branches which are
969   // not a TClonesArray
970
971
972   TString branchname;
973   TIter next(fESDObjects);
974   const Int_t kSplitlevel = 99; // default value in TTree::Branch()
975   const Int_t kBufsize = 32000; // default value in TTree::Branch()
976   TObject *obj = 0;
977
978   while ((obj = next())) {
979     branchname.Form("%s", obj->GetName());
980     if ((kSplitlevel > 1) &&  !obj->InheritsFrom(TClonesArray::Class())) {
981       branchname += ".";
982     }
983     tree->Bronch(branchname, obj->ClassName(), fESDObjects->GetObjectRef(obj),
984                  kBufsize, kSplitlevel - 1);
985   }
986
987 }
988
989
990 void AliESDEvent::ReadFromTree(TTree *tree){
991 //
992 // Connect the ESDEvent to a tree
993 //
994   if(!tree){
995     Printf("%s %d AliESDEvent::ReadFromTree() Zero Pointer to Tree \n",(char*)__FILE__,__LINE__);
996     return;
997   }
998   // load the TTree
999   if(!tree->GetTree())tree->LoadTree(0);
1000
1001   // if we find the "ESD" branch on the tree we do have the old structure
1002   if(tree->GetBranch("ESD")) {
1003     char ** address  = (char **)(tree->GetBranch("ESD")->GetAddress());
1004     // do we have the friend branch
1005     TBranch * esdFB = tree->GetBranch("ESDfriend.");
1006     char ** addressF = 0;
1007     if(esdFB)addressF = (char **)(esdFB->GetAddress());
1008     if (!address) {
1009       printf("%s %d AliESDEvent::ReadFromTree() Reading old Tree \n",(char*)__FILE__,__LINE__);
1010       tree->SetBranchAddress("ESD",       &fESDOld);
1011       if(esdFB){
1012         tree->SetBranchAddress("ESDfriend.",&fESDFriendOld);
1013       }
1014     } else {
1015       printf("%s %d AliESDEvent::ReadFromTree() Reading old Tree \n",(char*)__FILE__,__LINE__);
1016       printf("%s %d Branch already connected. Using existing branch address. \n",(char*)__FILE__,__LINE__);
1017       fESDOld       = (AliESD*)       (*address);
1018       // addressF can still be 0, since branch needs to switched on
1019       if(addressF)fESDFriendOld = (AliESDfriend*) (*addressF);
1020     }
1021                                        
1022     //  have already connected the old ESD structure... ?
1023     // reuse also the pointer of the AlliESDEvent
1024     // otherwise create new ones
1025     TList* connectedList = (TList*) (tree->GetUserInfo()->FindObject("ESDObjectsConnectedToTree"));
1026   
1027     if(connectedList){
1028       // If connected use the connected list of objects
1029       if(fESDObjects!= connectedList){
1030         // protect when called twice 
1031         fESDObjects->Delete();
1032         fESDObjects = connectedList;
1033       }
1034       GetStdContent(); 
1035
1036       
1037       // The pointer to the friend changes when called twice via InitIO
1038       // since AliESDEvent is deleted
1039       TObject* oldf = FindListObject("AliESDfriend");
1040       TObject* newf = 0;
1041       if(addressF){
1042         newf = (TObject*)*addressF;
1043       }
1044       if(newf!=0&&oldf!=newf){
1045         // remove the old reference
1046         // Should we also delete it? Or is this handled in TTree I/O
1047         // since it is created by the first SetBranchAddress
1048         fESDObjects->Remove(oldf);
1049         // add the new one 
1050         fESDObjects->Add(newf);
1051       }
1052       
1053       fConnected = true;
1054       return;
1055     }
1056     // else...    
1057     CreateStdContent(); // create for copy
1058     // if we have the esdfriend add it, so we always can access it via the userinfo
1059     if(fESDFriendOld)AddObject(fESDFriendOld);
1060     // we are not owner of the list objects 
1061     // must not delete it
1062     fESDObjects->SetOwner(kFALSE);
1063     fESDObjects->SetName("ESDObjectsConnectedToTree");
1064     tree->GetUserInfo()->Add(fESDObjects);
1065     fConnected = true;
1066     return;
1067   }
1068   
1069   delete fESDOld;
1070   fESDOld = 0;
1071   // Try to find AliESDEvent
1072   AliESDEvent *esdEvent = 0;
1073   esdEvent = (AliESDEvent*)tree->GetTree()->GetUserInfo()->FindObject("AliESDEvent");
1074   if(esdEvent){   
1075       // Check if already connected to tree
1076     TList* connectedList = (TList*) (tree->GetUserInfo()->FindObject("ESDObjectsConnectedToTree"));
1077     if (connectedList) {
1078       // If connected use the connected list if objects
1079       fESDObjects->Delete();
1080       fESDObjects = connectedList;
1081       GetStdContent(); 
1082       fConnected = true;
1083       return;
1084     }
1085     // Connect to tree
1086     // prevent a memory leak when reading back the TList
1087     delete fESDObjects;
1088     fESDObjects = 0;
1089     // create a new TList from the UserInfo TList... 
1090     // copy constructor does not work...
1091     fESDObjects = (TList*)(esdEvent->GetList()->Clone());
1092     fESDObjects->SetOwner(kFALSE);
1093     if(fESDObjects->GetEntries()<kESDListN){
1094       printf("%s %d AliESDEvent::ReadFromTree() TList contains less than the standard contents %d < %d \n",
1095              (char*)__FILE__,__LINE__,fESDObjects->GetEntries(),kESDListN);
1096     }
1097     // set the branch addresses
1098     TIter next(fESDObjects);
1099     TNamed *el;
1100     while((el=(TNamed*)next())){
1101       TString bname(el->GetName());
1102       if(bname.CompareTo("AliESDfriend")==0)
1103         {
1104           // AliESDfriend does not have a name ...
1105           tree->SetBranchAddress("ESDfriend.",fESDObjects->GetObjectRef(el));
1106         }
1107       else{
1108         // check if branch exists under this Name
1109         TBranch *br = tree->GetBranch(bname.Data());
1110         if(br){
1111           tree->SetBranchAddress(bname.Data(),fESDObjects->GetObjectRef(el));
1112         }
1113         else{
1114           br = tree->GetBranch(Form("%s.",bname.Data()));
1115           if(br){
1116             tree->SetBranchAddress(Form("%s.",bname.Data()),fESDObjects->GetObjectRef(el));
1117           }
1118           else{
1119             printf("%s %d AliESDEvent::ReadFromTree() No Branch found with Name %s or %s. \n",
1120                    (char*)__FILE__,__LINE__,bname.Data(),bname.Data());
1121           }
1122
1123         }
1124       }
1125     }
1126     GetStdContent();
1127     // when reading back we are not owner of the list 
1128     // must not delete it
1129     fESDObjects->SetOwner(kFALSE);
1130     fESDObjects->SetName("ESDObjectsConnectedToTree");
1131     // we are not owner of the list objects 
1132     // must not delete it
1133     tree->GetUserInfo()->Add(fESDObjects);
1134     fConnected = true;
1135   }// no esdEvent
1136   else {
1137     // we can't get the list from the user data, create standard content
1138     // and set it by hand (no ESDfriend at the moment
1139     CreateStdContent();
1140     TIter next(fESDObjects);
1141     TNamed *el;
1142     while((el=(TNamed*)next())){
1143       TString bname(el->GetName());    
1144       tree->SetBranchAddress(bname.Data(),fESDObjects->GetObjectRef(el));
1145     }
1146     GetStdContent();
1147     // when reading back we are not owner of the list 
1148     // must not delete it
1149     fESDObjects->SetOwner(kFALSE);
1150   }
1151 }
1152
1153
1154 void AliESDEvent::CopyFromOldESD()
1155 {
1156   // Method which copies over everthing from the old esd structure to the 
1157   // new  
1158   if(fESDOld){
1159     ResetStdContent();
1160      // Run
1161     SetRunNumber(fESDOld->GetRunNumber());
1162     SetPeriodNumber(fESDOld->GetPeriodNumber());
1163     SetMagneticField(fESDOld->GetMagneticField());
1164   
1165     // leave out diamond ...
1166     // SetDiamond(const AliESDVertex *vertex) { fESDRun->SetDiamond(vertex);}
1167
1168     // header
1169     SetTriggerMask(fESDOld->GetTriggerMask());
1170     SetOrbitNumber(fESDOld->GetOrbitNumber());
1171     SetTimeStamp(fESDOld->GetTimeStamp());
1172     SetEventType(fESDOld->GetEventType());
1173     SetEventNumberInFile(fESDOld->GetEventNumberInFile());
1174     SetBunchCrossNumber(fESDOld->GetBunchCrossNumber());
1175     SetTriggerCluster(fESDOld->GetTriggerCluster());
1176
1177     // ZDC
1178
1179     SetZDC(fESDOld->GetZDCN1Energy(),
1180            fESDOld->GetZDCP1Energy(),
1181            fESDOld->GetZDCEMEnergy(),
1182            0,
1183            fESDOld->GetZDCN2Energy(),
1184            fESDOld->GetZDCP2Energy(),
1185            fESDOld->GetZDCParticipants());
1186
1187     // FMD
1188     
1189     if(fESDOld->GetFMDData())SetFMDData(fESDOld->GetFMDData());
1190
1191     // T0
1192
1193     SetT0zVertex(fESDOld->GetT0zVertex());
1194     SetT0(fESDOld->GetT0());
1195     //  leave amps out
1196
1197     // VZERO
1198     if (fESDOld->GetVZEROData()) SetVZEROData(fESDOld->GetVZEROData());
1199
1200     if(fESDOld->GetVertex())SetPrimaryVertexSPD(fESDOld->GetVertex());
1201
1202     if(fESDOld->GetPrimaryVertex())SetPrimaryVertex(fESDOld->GetPrimaryVertex());
1203
1204     if(fESDOld->GetMultiplicity())SetMultiplicity(fESDOld->GetMultiplicity());
1205
1206     for(int i = 0;i<fESDOld->GetNumberOfTracks();i++){
1207       AddTrack(fESDOld->GetTrack(i));
1208     }
1209
1210     for(int i = 0;i<fESDOld->GetNumberOfMuonTracks();i++){
1211       AddMuonTrack(fESDOld->GetMuonTrack(i));
1212     }
1213
1214     for(int i = 0;i<fESDOld->GetNumberOfPmdTracks();i++){
1215       AddPmdTrack(fESDOld->GetPmdTrack(i));
1216     }
1217
1218     for(int i = 0;i<fESDOld->GetNumberOfTrdTracks();i++){
1219       AddTrdTrack(fESDOld->GetTrdTrack(i));
1220     }
1221
1222     for(int i = 0;i<fESDOld->GetNumberOfV0s();i++){
1223       AddV0(fESDOld->GetV0(i));
1224     }
1225
1226     for(int i = 0;i<fESDOld->GetNumberOfCascades();i++){
1227       AddCascade(fESDOld->GetCascade(i));
1228     }
1229
1230     for(int i = 0;i<fESDOld->GetNumberOfKinks();i++){
1231       AddKink(fESDOld->GetKink(i));
1232     }
1233
1234
1235     for(int i = 0;i<fESDOld->GetNumberOfCaloClusters();i++){
1236       AddCaloCluster(fESDOld->GetCaloCluster(i));
1237     }
1238
1239   }// if fesdold
1240 }
1241
1242
1243