1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
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.
23 AliESDEvent *ev= new AliESDEvent();
24 ev->ReadFromTree(esdTree);
26 for (Int_t i=0; i<nev; i++) {
28 if(ev->GetAliESDOld())ev->CopyFromOldESD();
30 // The AliESDInputHAndler does this automatically for you
32 // Origin: Christian Klein-Boesing, CERN, Christian.Klein-Boesing@cern.ch
33 //-----------------------------------------------------------------
36 #include "TRefArray.h"
39 #include "AliESDEvent.h"
40 #include "AliESDfriend.h"
41 #include "AliESDVZERO.h"
42 #include "AliESDFMD.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"
59 #include "AliESDFMD.h"
60 #include "AliESDVZERO.h"
61 #include "AliMultiplicity.h"
62 #include "AliRawDataErrorLog.h"
69 // here we define the names, some classes are no TNamed, therefore the classnames
71 const char* AliESDEvent::fgkESDListName[kESDListN] = {"AliESDRun",
92 "AliRawDataErrorLogs"};
93 //______________________________________________________________________________
94 AliESDEvent::AliESDEvent():
96 fESDObjects(new TList()),
116 fEMCALCells(0), fPHOSCells(0),
122 fFirstEMCALCluster(-1),
124 fFirstPHOSCluster(-1)
127 //______________________________________________________________________________
128 AliESDEvent::AliESDEvent(const AliESDEvent& esd):
130 fESDObjects(new TList()),
131 fESDRun(new AliESDRun(*esd.fESDRun)),
132 fHeader(new AliESDHeader(*esd.fHeader)),
133 fESDZDC(new AliESDZDC(*esd.fESDZDC)),
134 fESDFMD(new AliESDFMD(*esd.fESDFMD)),
135 fESDVZERO(new AliESDVZERO(*esd.fESDVZERO)),
136 fESDTZERO(new AliESDTZERO(*esd.fESDTZERO)),
137 fSPDVertex(new AliESDVertex(*esd.fSPDVertex)),
138 fPrimaryVertex(new AliESDVertex(*esd.fPrimaryVertex)),
139 fSPDMult(new AliMultiplicity(*esd.fSPDMult)),
140 fPHOSTrigger(new AliESDCaloTrigger(*esd.fPHOSTrigger)),
141 fEMCALTrigger(new AliESDCaloTrigger(*esd.fEMCALTrigger)),
142 fTracks(new TClonesArray(*esd.fTracks)),
143 fMuonTracks(new TClonesArray(*esd.fMuonTracks)),
144 fPmdTracks(new TClonesArray(*esd.fPmdTracks)),
145 fTrdTracks(new TClonesArray(*esd.fTrdTracks)),
146 fV0s(new TClonesArray(*esd.fV0s)),
147 fCascades(new TClonesArray(*esd.fCascades)),
148 fKinks(new TClonesArray(*esd.fKinks)),
149 fCaloClusters(new TClonesArray(*esd.fCaloClusters)),
150 fEMCALCells(new AliESDCaloCells(*esd.fEMCALCells)),
151 fPHOSCells(new AliESDCaloCells(*esd.fPHOSCells)),
152 fErrorLogs(new TClonesArray(*esd.fErrorLogs)),
153 fESDOld(new AliESD(*esd.fESDOld)),
154 fESDFriendOld(new AliESDfriend(*esd.fESDFriendOld)),
155 fConnected(esd.fConnected),
156 fEMCALClusters(esd.fEMCALClusters),
157 fFirstEMCALCluster(esd.fFirstEMCALCluster),
158 fPHOSClusters(esd.fPHOSClusters),
159 fFirstPHOSCluster(esd.fFirstPHOSCluster)
162 // CKB init in the constructor list and only add here ...
167 AddObject(fESDVZERO);
168 AddObject(fESDTZERO);
169 AddObject(fSPDVertex);
170 AddObject(fPrimaryVertex);
172 AddObject(fPHOSTrigger);
173 AddObject(fEMCALTrigger);
175 AddObject(fMuonTracks);
176 AddObject(fPmdTracks);
177 AddObject(fTrdTracks);
179 AddObject(fCascades);
181 AddObject(fCaloClusters);
182 AddObject(fEMCALCells);
183 AddObject(fPHOSCells);
184 AddObject(fErrorLogs);
190 //______________________________________________________________________________
191 AliESDEvent & AliESDEvent::operator=(const AliESDEvent& source) {
193 // Assignment operator
195 if(&source == this) return *this;
196 AliVEvent::operator=(source);
198 fESDRun = new AliESDRun(*source.fESDRun);
199 fHeader = new AliESDHeader(*source.fHeader);
200 fESDZDC = new AliESDZDC(*source.fESDZDC);
201 fESDFMD = new AliESDFMD(*source.fESDFMD);
202 fESDVZERO = new AliESDVZERO(*source.fESDVZERO);
203 fESDTZERO = new AliESDTZERO(*source.fESDTZERO);
204 fSPDVertex = new AliESDVertex(*source.fSPDVertex);
205 fPrimaryVertex = new AliESDVertex(*source.fPrimaryVertex);
206 fSPDMult = new AliMultiplicity(*source.fSPDMult);
207 fPHOSTrigger = new AliESDCaloTrigger(*source.fPHOSTrigger);
208 fEMCALTrigger = new AliESDCaloTrigger(*source.fEMCALTrigger);
209 fTracks = new TClonesArray(*source.fTracks);
210 fMuonTracks = new TClonesArray(*source.fMuonTracks);
211 fPmdTracks = new TClonesArray(*source.fPmdTracks);
212 fTrdTracks = new TClonesArray(*source.fTrdTracks);
213 fV0s = new TClonesArray(*source.fV0s);
214 fCascades = new TClonesArray(*source.fCascades);
215 fKinks = new TClonesArray(*source.fKinks);
216 fCaloClusters = new TClonesArray(*source.fCaloClusters);
217 fEMCALCells = new AliESDCaloCells(*source.fEMCALCells);
218 fPHOSCells = new AliESDCaloCells(*source.fPHOSCells);
219 fErrorLogs = new TClonesArray(*source.fErrorLogs);
220 fESDOld = new AliESD(*source.fESDOld);
221 fESDFriendOld = new AliESDfriend(*source.fESDFriendOld);
223 // or AddObject( fESDZDC = new AliESDZDC(*source.fESDZDC));
225 fESDObjects = new TList();
230 AddObject(fESDVZERO);
231 AddObject(fESDTZERO);
232 AddObject(fSPDVertex);
233 AddObject(fPrimaryVertex);
235 AddObject(fPHOSTrigger);
236 AddObject(fEMCALTrigger);
238 AddObject(fMuonTracks);
239 AddObject(fPmdTracks);
240 AddObject(fTrdTracks);
242 AddObject(fCascades);
244 AddObject(fCaloClusters);
245 AddObject(fEMCALCells);
246 AddObject(fPHOSCells);
247 AddObject(fErrorLogs);
249 fConnected = source.fConnected;
250 fEMCALClusters = source.fEMCALClusters;
251 fFirstEMCALCluster = source.fFirstEMCALCluster;
252 fPHOSClusters = source.fPHOSClusters;
253 fFirstPHOSCluster = source.fFirstPHOSCluster;
262 //______________________________________________________________________________
263 AliESDEvent::~AliESDEvent()
266 // Standard destructor
269 // everthing on the list gets deleted automatically
272 if(fESDObjects&&!fConnected)
281 //______________________________________________________________________________
282 void AliESDEvent::Reset()
286 // Reset the standard contents
288 if(fESDOld)fESDOld->Reset();
289 // reset for the friends...
291 fESDFriendOld->~AliESDfriend();
292 new (fESDFriendOld) AliESDfriend();
294 // for new data we have to fetch the Pointer from the list
295 AliESDfriend *fr = (AliESDfriend*)FindListObject("AliESDfriend");
297 // delete the content
299 // make a new valid ESDfriend at the same place
300 new (fr) AliESDfriend();
303 // call reset for user supplied data?
306 void AliESDEvent::ResetStdContent()
308 // Reset the standard contents
309 if(fESDRun) fESDRun->Reset();
310 if(fHeader) fHeader->Reset();
311 if(fESDZDC) fESDZDC->Reset();
312 if(fESDFMD) fESDFMD->Clear(); // why clear.... need consistend names
314 // reset by callin d'to /c'tor keep the pointer
315 fESDVZERO->~AliESDVZERO();
316 new (fESDVZERO) AliESDVZERO();
318 if(fESDTZERO) fESDTZERO->Reset();
319 // CKB no clear/reset implemented
321 fSPDVertex->~AliESDVertex();
322 new (fSPDVertex) AliESDVertex();
323 fSPDVertex->SetName(fgkESDListName[kSPDVertex]);
326 fPrimaryVertex->~AliESDVertex();
327 new (fPrimaryVertex) AliESDVertex();
328 fPrimaryVertex->SetName(fgkESDListName[kPrimaryVertex]);
331 fSPDMult->~AliMultiplicity();
332 new (fSPDMult) AliMultiplicity();
334 if(fPHOSTrigger)fPHOSTrigger->Reset();
335 if(fEMCALTrigger)fEMCALTrigger->Reset();
336 if(fTracks)fTracks->Delete();
337 if(fMuonTracks)fMuonTracks->Delete();
338 if(fPmdTracks)fPmdTracks->Delete();
339 if(fTrdTracks)fTrdTracks->Delete();
340 if(fV0s)fV0s->Delete();
341 if(fCascades)fCascades->Delete();
342 if(fKinks)fKinks->Delete();
343 if(fCaloClusters)fCaloClusters->Delete();
344 if(fPHOSCells)fPHOSCells->DeleteContainer();
345 if(fEMCALCells)fEMCALCells->DeleteContainer();
346 if(fErrorLogs) fErrorLogs->Delete();
348 // don't reset fconnected fConnected and the list
351 fFirstEMCALCluster=-1;
353 fFirstPHOSCluster=-1;
357 Int_t AliESDEvent::AddV0(const AliESDv0 *v) {
361 TClonesArray &fv = *fV0s;
362 Int_t idx=fV0s->GetEntriesFast();
363 new(fv[idx]) AliESDv0(*v);
367 //______________________________________________________________________________
368 void AliESDEvent::Print(Option_t *) const
371 // Print header information of the event
373 printf("ESD run information\n");
374 printf("Event # in file %d Bunch crossing # %d Orbit # %d Period # %d Run # %d Trigger %lld Magnetic field %f \n",
375 GetEventNumberInFile(),
376 GetBunchCrossNumber(),
381 GetMagneticField() );
382 printf("Vertex: (%.4f +- %.4f, %.4f +- %.4f, %.4f +- %.4f) cm\n",
383 fPrimaryVertex->GetXv(), fPrimaryVertex->GetXRes(),
384 fPrimaryVertex->GetYv(), fPrimaryVertex->GetYRes(),
385 fPrimaryVertex->GetZv(), fPrimaryVertex->GetZRes());
386 printf("Mean vertex in RUN: X=%.4f Y=%.4f cm\n",
387 GetDiamondX(),GetDiamondY());
388 printf("SPD Multiplicity. Number of tracklets %d \n",
389 fSPDMult->GetNumberOfTracklets());
390 printf("Number of tracks: \n");
391 printf(" charged %d\n", GetNumberOfTracks());
392 printf(" muon %d\n", GetNumberOfMuonTracks());
393 printf(" pmd %d\n", GetNumberOfPmdTracks());
394 printf(" trd %d\n", GetNumberOfTrdTracks());
395 printf(" v0 %d\n", GetNumberOfV0s());
396 printf(" cascades %d\n", GetNumberOfCascades());
397 printf(" kinks %d\n", GetNumberOfKinks());
398 if(fPHOSCells)printf(" PHOSCells %d\n", fPHOSCells->GetNumberOfCells());
399 else printf(" PHOSCells not in the Event\n");
400 if(fEMCALCells)printf(" EMCALCells %d\n", fEMCALCells->GetNumberOfCells());
401 else printf(" EMCALCells not in the Event\n");
402 printf(" CaloClusters %d\n", GetNumberOfCaloClusters());
403 printf(" phos %d\n", GetNumberOfPHOSClusters());
404 printf(" emcal %d\n", GetNumberOfEMCALClusters());
405 printf(" FMD %s\n", (fESDFMD ? "yes" : "no"));
406 printf(" VZERO %s\n", (fESDVZERO ? "yes" : "no"));
411 void AliESDEvent::SetESDfriend(const AliESDfriend *ev) {
413 // Attaches the complementary info to the ESD
417 // to be sure that we set the tracks also
418 // in case of old esds
419 // if(fESDOld)CopyFromOldESD();
421 Int_t ntrk=ev->GetNumberOfTracks();
423 for (Int_t i=0; i<ntrk; i++) {
424 const AliESDfriendTrack *f=ev->GetTrack(i);
425 GetTrack(i)->SetFriendTrack(f);
429 Bool_t AliESDEvent::RemoveKink(Int_t rm) {
430 // ---------------------------------------------------------
431 // Remove a kink candidate and references to it from ESD,
432 // if this candidate does not come from a reconstructed decay
433 // Not yet implemented...
434 // ---------------------------------------------------------
435 Int_t last=GetNumberOfKinks()-1;
436 if ((rm<0)||(rm>last)) return kFALSE;
441 Bool_t AliESDEvent::RemoveV0(Int_t rm) {
442 // ---------------------------------------------------------
443 // Remove a V0 candidate and references to it from ESD,
444 // if this candidate does not come from a reconstructed decay
445 // ---------------------------------------------------------
446 Int_t last=GetNumberOfV0s()-1;
447 if ((rm<0)||(rm>last)) return kFALSE;
449 AliESDv0 *v0=GetV0(rm);
450 Int_t idxP=v0->GetPindex(), idxN=v0->GetNindex();
453 Int_t lastIdxP=v0->GetPindex(), lastIdxN=v0->GetNindex();
457 // Check if this V0 comes from a reconstructed decay
458 Int_t ncs=GetNumberOfCascades();
459 for (Int_t n=0; n<ncs; n++) {
460 AliESDcascade *cs=GetCascade(n);
462 Int_t csIdxP=cs->GetPindex();
463 Int_t csIdxN=cs->GetNindex();
466 if (idxN==csIdxN) return kFALSE;
468 if (csIdxP==lastIdxP)
469 if (csIdxN==lastIdxN) used++;
472 //Replace the removed V0 with the last V0
473 TClonesArray &a=*fV0s;
474 delete a.RemoveAt(rm);
476 if (rm==last) return kTRUE;
478 //v0 is pointing to the last V0 candidate...
479 new (a[rm]) AliESDv0(*v0);
480 delete a.RemoveAt(last);
482 if (!used) return kTRUE;
485 // Remap the indices of the daughters of reconstructed decays
486 for (Int_t n=0; n<ncs; n++) {
487 AliESDcascade *cs=GetCascade(n);
490 Int_t csIdxP=cs->GetPindex();
491 Int_t csIdxN=cs->GetNindex();
493 if (csIdxP==lastIdxP)
494 if (csIdxN==lastIdxN) {
495 cs->AliESDv0::SetIndex(1,idxP);
496 cs->AliESDv0::SetIndex(0,idxN);
498 if (!used) return kTRUE;
505 Bool_t AliESDEvent::RemoveTrack(Int_t rm) {
506 // ---------------------------------------------------------
507 // Remove a track and references to it from ESD,
508 // if this track does not come from a reconstructed decay
509 // ---------------------------------------------------------
510 Int_t last=GetNumberOfTracks()-1;
511 if ((rm<0)||(rm>last)) return kFALSE;
515 // Check if this track comes from a reconstructed decay
516 Int_t nv0=GetNumberOfV0s();
517 for (Int_t n=0; n<nv0; n++) {
518 AliESDv0 *v0=GetV0(n);
520 Int_t idx=v0->GetNindex();
521 if (rm==idx) return kFALSE;
522 if (idx==last) used++;
525 if (rm==idx) return kFALSE;
526 if (idx==last) used++;
529 Int_t ncs=GetNumberOfCascades();
530 for (Int_t n=0; n<ncs; n++) {
531 AliESDcascade *cs=GetCascade(n);
533 Int_t idx=cs->GetIndex();
534 if (rm==idx) return kFALSE;
535 if (idx==last) used++;
538 Int_t nkn=GetNumberOfKinks();
539 for (Int_t n=0; n<nkn; n++) {
540 AliESDkink *kn=GetKink(n);
542 Int_t idx=kn->GetIndex(0);
543 if (rm==idx) return kFALSE;
544 if (idx==last) used++;
547 if (rm==idx) return kFALSE;
548 if (idx==last) used++;
552 //Replace the removed track with the last track
553 TClonesArray &a=*fTracks;
554 delete a.RemoveAt(rm);
556 if (rm==last) return kTRUE;
558 AliESDtrack *t=GetTrack(last);
560 new (a[rm]) AliESDtrack(*t);
561 delete a.RemoveAt(last);
563 if (!used) return kTRUE;
566 // Remap the indices of the daughters of reconstructed decays
567 for (Int_t n=0; n<nv0; n++) {
568 AliESDv0 *v0=GetV0(n);
569 if (v0->GetIndex(0)==last) {
572 if (!used) return kTRUE;
574 if (v0->GetIndex(1)==last) {
577 if (!used) return kTRUE;
581 for (Int_t n=0; n<ncs; n++) {
582 AliESDcascade *cs=GetCascade(n);
583 if (cs->GetIndex()==last) {
586 if (!used) return kTRUE;
590 for (Int_t n=0; n<nkn; n++) {
591 AliESDkink *kn=GetKink(n);
592 if (kn->GetIndex(0)==last) {
595 if (!used) return kTRUE;
597 if (kn->GetIndex(1)==last) {
600 if (!used) return kTRUE;
608 Bool_t AliESDEvent::Clean(Float_t *cleanPars) {
610 // Remove the data which are not needed for the physics analysis.
612 // 1) Cleaning the V0 candidates
613 // ---------------------------
614 // If the cosine of the V0 pointing angle "csp" and
615 // the DCA between the daughter tracks "dca" does not satisfy
618 // csp > cleanPars[1] + dca/cleanPars[0]*(1.- cleanPars[1])
620 // an attempt to remove this V0 candidate from ESD is made.
622 // The V0 candidate gets removed if it does not belong to any
623 // recosntructed cascade decay
625 // 12.11.2007, optimal values: cleanPars[0]=0.5, cleanPars[1]=0.999
627 // 2) Cleaning the tracks
628 // ----------------------
629 // If track's transverse parameter is larger than cleanPars[2]
631 // track's longitudinal parameter is larger than cleanPars[3]
632 // an attempt to remove this track from ESD is made.
634 // The track gets removed if it does not come
635 // from a reconstructed decay
639 Float_t dcaMax=cleanPars[0];
640 Float_t cspMin=cleanPars[1];
642 Int_t nV0s=GetNumberOfV0s();
643 for (Int_t i=nV0s-1; i>=0; i--) {
644 AliESDv0 *v0=GetV0(i);
646 Float_t dca=v0->GetDcaV0Daughters();
647 Float_t csp=v0->GetV0CosineOfPointingAngle();
648 Float_t cspcut=cspMin + dca/dcaMax*(1.-cspMin);
650 if (RemoveV0(i)) rc=kTRUE;
655 Float_t dmax=cleanPars[2], zmax=cleanPars[3];
657 const AliESDVertex *vertex=GetVertex();
658 Bool_t vtxOK=vertex->GetStatus();
660 Int_t nTracks=GetNumberOfTracks();
661 for (Int_t i=nTracks-1; i>=0; i--) {
662 AliESDtrack *track=GetTrack(i);
663 Float_t xy,z; track->GetImpactParameters(xy,z);
664 if ((TMath::Abs(xy) > dmax) || (vtxOK && (TMath::Abs(z) > zmax))) {
665 if (RemoveTrack(i)) rc=kTRUE;
672 Int_t AliESDEvent::AddTrack(const AliESDtrack *t)
675 TClonesArray &ftr = *fTracks;
676 AliESDtrack * track = new(ftr[fTracks->GetEntriesFast()])AliESDtrack(*t);
677 track->SetID(fTracks->GetEntriesFast()-1);
678 return track->GetID();
681 void AliESDEvent::AddMuonTrack(const AliESDMuonTrack *t)
683 TClonesArray &fmu = *fMuonTracks;
684 new(fmu[fMuonTracks->GetEntriesFast()]) AliESDMuonTrack(*t);
687 void AliESDEvent::AddPmdTrack(const AliESDPmdTrack *t)
689 TClonesArray &fpmd = *fPmdTracks;
690 new(fpmd[fPmdTracks->GetEntriesFast()]) AliESDPmdTrack(*t);
693 void AliESDEvent::AddTrdTrack(const AliESDTrdTrack *t)
695 TClonesArray &ftrd = *fTrdTracks;
696 new(ftrd[fTrdTracks->GetEntriesFast()]) AliESDTrdTrack(*t);
702 Int_t AliESDEvent::AddKink(const AliESDkink *c)
705 TClonesArray &fk = *fKinks;
706 AliESDkink * kink = new(fk[fKinks->GetEntriesFast()]) AliESDkink(*c);
707 kink->SetID(fKinks->GetEntriesFast()); // CKB different from the other imps..
708 return fKinks->GetEntriesFast()-1;
712 void AliESDEvent::AddCascade(const AliESDcascade *c)
714 TClonesArray &fc = *fCascades;
715 new(fc[fCascades->GetEntriesFast()]) AliESDcascade(*c);
719 Int_t AliESDEvent::AddCaloCluster(const AliESDCaloCluster *c)
722 TClonesArray &fc = *fCaloClusters;
723 AliESDCaloCluster *clus = new(fc[fCaloClusters->GetEntriesFast()]) AliESDCaloCluster(*c);
724 clus->SetID(fCaloClusters->GetEntriesFast()-1);
725 return fCaloClusters->GetEntriesFast()-1;
729 void AliESDEvent::AddRawDataErrorLog(const AliRawDataErrorLog *log) {
730 TClonesArray &errlogs = *fErrorLogs;
731 new(errlogs[errlogs.GetEntriesFast()]) AliRawDataErrorLog(*log);
734 void AliESDEvent::SetVertex(const AliESDVertex *vertex)
736 // Set the SPD vertex
737 // use already allocated space
739 *fSPDVertex = *vertex;
740 fSPDVertex->SetName(fgkESDListName[kSPDVertex]);
744 void AliESDEvent::SetPrimaryVertex(const AliESDVertex *vertex)
746 // Set the primary vertex
747 // use already allocated space
749 *fPrimaryVertex = *vertex;
750 fPrimaryVertex->SetName(fgkESDListName[kPrimaryVertex]);
754 void AliESDEvent::SetMultiplicity(const AliMultiplicity *mul)
756 // Set the SPD Multiplicity
763 void AliESDEvent::SetFMDData(AliESDFMD * obj)
765 // use already allocated space
771 void AliESDEvent::SetVZEROData(AliESDVZERO * obj)
773 // use already allocated space
778 void AliESDEvent::GetESDfriend(AliESDfriend *ev) const
781 // Extracts the complementary info from the ESD
785 Int_t ntrk=GetNumberOfTracks();
787 for (Int_t i=0; i<ntrk; i++) {
788 AliESDtrack *t=GetTrack(i);
789 const AliESDfriendTrack *f=t->GetFriendTrack();
792 t->ReleaseESDfriendTrack();// Not to have two copies of "friendTrack"
798 void AliESDEvent::AddObject(TObject* obj)
800 // Add an object to the list of object.
801 // Please be aware that in order to increase performance you should
802 // refrain from using TObjArrays (if possible). Use TClonesArrays, instead.
803 fESDObjects->SetOwner(kTRUE);
804 fESDObjects->AddLast(obj);
808 void AliESDEvent::GetStdContent()
810 // set pointers for standard content
811 // get by name much safer and not a big overhead since not called very often
813 fESDRun = (AliESDRun*)fESDObjects->FindObject(fgkESDListName[kESDRun]);
814 fHeader = (AliESDHeader*)fESDObjects->FindObject(fgkESDListName[kHeader]);
815 fESDZDC = (AliESDZDC*)fESDObjects->FindObject(fgkESDListName[kESDZDC]);
816 fESDFMD = (AliESDFMD*)fESDObjects->FindObject(fgkESDListName[kESDFMD]);
817 fESDVZERO = (AliESDVZERO*)fESDObjects->FindObject(fgkESDListName[kESDVZERO]);
818 fESDTZERO = (AliESDTZERO*)fESDObjects->FindObject(fgkESDListName[kESDTZERO]);
819 fSPDVertex = (AliESDVertex*)fESDObjects->FindObject(fgkESDListName[kSPDVertex]);
820 fPrimaryVertex = (AliESDVertex*)fESDObjects->FindObject(fgkESDListName[kPrimaryVertex]);
821 fSPDMult = (AliMultiplicity*)fESDObjects->FindObject(fgkESDListName[kSPDMult]);
822 fPHOSTrigger = (AliESDCaloTrigger*)fESDObjects->FindObject(fgkESDListName[kPHOSTrigger]);
823 fEMCALTrigger = (AliESDCaloTrigger*)fESDObjects->FindObject(fgkESDListName[kEMCALTrigger]);
824 fTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kTracks]);
825 fMuonTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kMuonTracks]);
826 fPmdTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kPmdTracks]);
827 fTrdTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kTrdTracks]);
828 fV0s = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kV0s]);
829 fCascades = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kCascades]);
830 fKinks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kKinks]);
831 fCaloClusters = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kCaloClusters]);
832 fEMCALCells = (AliESDCaloCells*)fESDObjects->FindObject(fgkESDListName[kEMCALCells]);
833 fPHOSCells = (AliESDCaloCells*)fESDObjects->FindObject(fgkESDListName[kPHOSCells]);
834 fErrorLogs = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kErrorLogs]);
838 void AliESDEvent::SetStdNames(){
839 // Set the names of the standard contents
841 if(fESDObjects->GetEntries()==kESDListN){
842 for(int i = 0;i < fESDObjects->GetEntries();i++){
843 TObject *fObj = fESDObjects->At(i);
844 if(fObj->InheritsFrom("TNamed")){
845 ((TNamed*)fObj)->SetName(fgkESDListName[i]);
847 else if(fObj->InheritsFrom("TClonesArray")){
848 ((TClonesArray*)fObj)->SetName(fgkESDListName[i]);
853 printf("%s:%d SetStdNames() Wrong number of Std Entries \n",(char*)__FILE__,__LINE__);
857 void AliESDEvent::CreateStdContent()
859 // create the standard AOD content and set pointers
861 // create standard objects and add them to the TList of objects
862 AddObject(new AliESDRun());
863 AddObject(new AliESDHeader());
864 AddObject(new AliESDZDC());
865 AddObject(new AliESDFMD());
866 AddObject(new AliESDVZERO());
867 AddObject(new AliESDTZERO());
868 AddObject(new AliESDVertex());
869 AddObject(new AliESDVertex());
870 AddObject(new AliMultiplicity());
871 AddObject(new AliESDCaloTrigger());
872 AddObject(new AliESDCaloTrigger());
873 AddObject(new TClonesArray("AliESDtrack",0));
874 AddObject(new TClonesArray("AliESDMuonTrack",0));
875 AddObject(new TClonesArray("AliESDPmdTrack",0));
876 AddObject(new TClonesArray("AliESDTrdTrack",0));
877 AddObject(new TClonesArray("AliESDv0",0));
878 AddObject(new TClonesArray("AliESDcascade",0));
879 AddObject(new TClonesArray("AliESDkink",0));
880 AddObject(new TClonesArray("AliESDCaloCluster",0));
881 AddObject(new AliESDCaloCells());
882 AddObject(new AliESDCaloCells());
883 AddObject(new TClonesArray("AliRawDataErrorLog",0));
885 // check the order of the indices against enum...
889 // read back pointers
893 TObject* AliESDEvent::FindListObject(const char *name){
895 return fESDObjects->FindObject(name);
900 Int_t AliESDEvent::GetPHOSClusters(TRefArray *clusters) const
902 // fills the provided TRefArray with all found phos clusters
906 AliESDCaloCluster *cl = 0;
907 for (Int_t i = 0; i < GetNumberOfCaloClusters(); i++) {
909 if ( (cl = GetCaloCluster(i)) ) {
912 AliDebug(1,Form("IsPHOS cluster %d Size: %d \n",i,clusters->GetEntriesFast()));
916 return clusters->GetEntriesFast();
919 Int_t AliESDEvent::GetEMCALClusters(TRefArray *clusters) const
921 // fills the provided TRefArray with all found emcal clusters
925 AliESDCaloCluster *cl = 0;
926 for (Int_t i = 0; i < GetNumberOfCaloClusters(); i++) {
928 if ( (cl = GetCaloCluster(i)) ) {
931 AliDebug(1,Form("IsEMCAL cluster %d Size: %d \n",i,clusters->GetEntriesFast()));
935 return clusters->GetEntriesFast();
938 const void AliESDEvent::WriteToTree(TTree* tree) const {
939 // Book the branches as in TTree::Branch(TCollection*)
940 // but add a "." at the end of top level branches which are
941 // not a TClonesArray
945 TIter next(fESDObjects);
946 const Int_t splitlevel = 99; // default value in TTree::Branch()
947 const Int_t bufsize = 32000; // default value in TTree::Branch()
950 while ((obj = next())) {
951 branchname.Form("%s", obj->GetName());
952 if ((splitlevel > 1) && !obj->InheritsFrom(TClonesArray::Class())) {
955 tree->Bronch(branchname, obj->ClassName(), fESDObjects->GetObjectRef(obj),
956 bufsize, splitlevel - 1);
962 void AliESDEvent::ReadFromTree(TTree *tree){
965 Printf("%s %d AliESDEvent::ReadFromTree() Zero Pointer to Tree \n",(char*)__FILE__,__LINE__);
969 if(!tree->GetTree())tree->LoadTree(0);
971 // if we find the "ESD" branch on the tree we do have the old structure
972 if(tree->GetBranch("ESD")) {
973 char ** address = (char **)(tree->GetBranch("ESD")->GetAddress());
974 // do we have the friend branch
975 TBranch * esdFB = tree->GetBranch("ESDfriend.");
976 char ** addressF = 0;
977 if(esdFB)addressF = (char **)(esdFB->GetAddress());
979 printf("%s %d AliESDEvent::ReadFromTree() Reading old Tree \n",(char*)__FILE__,__LINE__);
980 tree->SetBranchAddress("ESD", &fESDOld);
982 tree->SetBranchAddress("ESDfriend.",&fESDFriendOld);
985 printf("%s %d AliESDEvent::ReadFromTree() Reading old Tree \n",(char*)__FILE__,__LINE__);
986 printf("%s %d Branch already connected. Using existing branch address. \n",(char*)__FILE__,__LINE__);
987 fESDOld = (AliESD*) (*address);
988 // addressF can still be 0, since branch needs to switched on
989 if(addressF)fESDFriendOld = (AliESDfriend*) (*addressF);
992 // have already connected the old ESD structure... ?
993 // reuse also the pointer of the AlliESDEvent
994 // otherwise create new ones
995 TList* connectedList = (TList*) (tree->GetUserInfo()->FindObject("ESDObjectsConnectedToTree"));
998 // If connected use the connected list of objects
999 if(fESDObjects!= connectedList){
1000 // protect when called twice
1001 fESDObjects->Delete();
1002 fESDObjects = connectedList;
1007 // The pointer to the friend changes when called twice via InitIO
1008 // since AliESDEvent is deleted
1009 TObject* oldf = FindListObject("AliESDfriend");
1012 newf = (TObject*)*addressF;
1014 if(newf!=0&&oldf!=newf){
1015 // remove the old reference
1016 // Should we also delete it? Or is this handled in TTree I/O
1017 // since it is created by the first SetBranchAddress
1018 fESDObjects->Remove(oldf);
1020 fESDObjects->Add(newf);
1027 CreateStdContent(); // create for copy
1028 // if we have the esdfriend add it, so we always can access it via the userinfo
1029 if(fESDFriendOld)AddObject(fESDFriendOld);
1030 // we are not owner of the list objects
1031 // must not delete it
1032 fESDObjects->SetOwner(kFALSE);
1033 fESDObjects->SetName("ESDObjectsConnectedToTree");
1034 tree->GetUserInfo()->Add(fESDObjects);
1041 // Try to find AliESDEvent
1042 AliESDEvent *esdEvent = 0;
1043 esdEvent = (AliESDEvent*)tree->GetTree()->GetUserInfo()->FindObject("AliESDEvent");
1045 // Check if already connected to tree
1046 TList* connectedList = (TList*) (tree->GetUserInfo()->FindObject("ESDObjectsConnectedToTree"));
1047 if (connectedList) {
1048 // If connected use the connected list if objects
1049 fESDObjects->Delete();
1050 fESDObjects = connectedList;
1056 // prevent a memory leak when reading back the TList
1059 // create a new TList from the UserInfo TList...
1060 // copy constructor does not work...
1061 fESDObjects = (TList*)(esdEvent->GetList()->Clone());
1062 fESDObjects->SetOwner(kFALSE);
1063 if(fESDObjects->GetEntries()<kESDListN){
1064 printf("%s %d AliESDEvent::ReadFromTree() TList contains less than the standard contents %d < %d \n",
1065 (char*)__FILE__,__LINE__,fESDObjects->GetEntries(),kESDListN);
1067 // set the branch addresses
1068 TIter next(fESDObjects);
1070 while((el=(TNamed*)next())){
1071 TString bname(el->GetName());
1072 if(bname.CompareTo("AliESDfriend")==0)
1074 // AliESDfriend does not have a name ...
1075 tree->SetBranchAddress("ESDfriend.",fESDObjects->GetObjectRef(el));
1078 // check if branch exists under this Name
1079 TBranch *br = tree->GetBranch(bname.Data());
1081 tree->SetBranchAddress(bname.Data(),fESDObjects->GetObjectRef(el));
1084 br = tree->GetBranch(Form("%s.",bname.Data()));
1086 tree->SetBranchAddress(Form("%s.",bname.Data()),fESDObjects->GetObjectRef(el));
1089 printf("%s %d AliESDEvent::ReadFromTree() No Branch found with Name %s or %s. \n",
1090 (char*)__FILE__,__LINE__,bname.Data(),bname.Data());
1097 // when reading back we are not owner of the list
1098 // must not delete it
1099 fESDObjects->SetOwner(kFALSE);
1100 fESDObjects->SetName("ESDObjectsConnectedToTree");
1101 // we are not owner of the list objects
1102 // must not delete it
1103 tree->GetUserInfo()->Add(fESDObjects);
1107 // we can't get the list from the user data, create standard content
1108 // and set it by hand (no ESDfriend at the moment
1110 TIter next(fESDObjects);
1112 while((el=(TNamed*)next())){
1113 TString bname(el->GetName());
1114 tree->SetBranchAddress(bname.Data(),fESDObjects->GetObjectRef(el));
1117 // when reading back we are not owner of the list
1118 // must not delete it
1119 fESDObjects->SetOwner(kFALSE);
1124 void AliESDEvent::CopyFromOldESD()
1126 // Method which copies over everthing from the old esd structure to the
1131 SetRunNumber(fESDOld->GetRunNumber());
1132 SetPeriodNumber(fESDOld->GetPeriodNumber());
1133 SetMagneticField(fESDOld->GetMagneticField());
1135 // leave out diamond ...
1136 // SetDiamond(const AliESDVertex *vertex) { fESDRun->SetDiamond(vertex);}
1139 SetTriggerMask(fESDOld->GetTriggerMask());
1140 SetOrbitNumber(fESDOld->GetOrbitNumber());
1141 SetTimeStamp(fESDOld->GetTimeStamp());
1142 SetEventType(fESDOld->GetEventType());
1143 SetEventNumberInFile(fESDOld->GetEventNumberInFile());
1144 SetBunchCrossNumber(fESDOld->GetBunchCrossNumber());
1145 SetTriggerCluster(fESDOld->GetTriggerCluster());
1149 SetZDC(fESDOld->GetZDCN1Energy(),
1150 fESDOld->GetZDCP1Energy(),
1151 fESDOld->GetZDCEMEnergy(),
1153 fESDOld->GetZDCN2Energy(),
1154 fESDOld->GetZDCP2Energy(),
1155 fESDOld->GetZDCParticipants());
1159 if(fESDOld->GetFMDData())SetFMDData(fESDOld->GetFMDData());
1163 SetT0zVertex(fESDOld->GetT0zVertex());
1164 SetT0(fESDOld->GetT0());
1168 if (fESDOld->GetVZEROData()) SetVZEROData(fESDOld->GetVZEROData());
1170 if(fESDOld->GetVertex())SetVertex(fESDOld->GetVertex());
1172 if(fESDOld->GetPrimaryVertex())SetPrimaryVertex(fESDOld->GetPrimaryVertex());
1174 if(fESDOld->GetMultiplicity())SetMultiplicity(fESDOld->GetMultiplicity());
1176 for(int i = 0;i<fESDOld->GetNumberOfTracks();i++){
1177 AddTrack(fESDOld->GetTrack(i));
1180 for(int i = 0;i<fESDOld->GetNumberOfMuonTracks();i++){
1181 AddMuonTrack(fESDOld->GetMuonTrack(i));
1184 for(int i = 0;i<fESDOld->GetNumberOfPmdTracks();i++){
1185 AddPmdTrack(fESDOld->GetPmdTrack(i));
1188 for(int i = 0;i<fESDOld->GetNumberOfTrdTracks();i++){
1189 AddTrdTrack(fESDOld->GetTrdTrack(i));
1192 for(int i = 0;i<fESDOld->GetNumberOfV0s();i++){
1193 AddV0(fESDOld->GetV0(i));
1196 for(int i = 0;i<fESDOld->GetNumberOfCascades();i++){
1197 AddCascade(fESDOld->GetCascade(i));
1200 for(int i = 0;i<fESDOld->GetNumberOfKinks();i++){
1201 AddKink(fESDOld->GetKink(i));
1205 for(int i = 0;i<fESDOld->GetNumberOfCaloClusters();i++){
1206 AddCaloCluster(fESDOld->GetCaloCluster(i));