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"
64 #include "AliESDACORDE.h"
69 // here we define the names, some classes are no TNamed, therefore the classnames
71 const char* AliESDEvent::fgkESDListName[kESDListN] = {"AliESDRun",
93 "AliRawDataErrorLogs",
96 //______________________________________________________________________________
97 AliESDEvent::AliESDEvent():
99 fESDObjects(new TList()),
121 fEMCALCells(0), fPHOSCells(0),
128 fFirstEMCALCluster(-1),
130 fFirstPHOSCluster(-1)
133 //______________________________________________________________________________
134 AliESDEvent::AliESDEvent(const AliESDEvent& esd):
136 fESDObjects(new TList()),
137 fESDRun(new AliESDRun(*esd.fESDRun)),
138 fHeader(new AliESDHeader(*esd.fHeader)),
139 fESDZDC(new AliESDZDC(*esd.fESDZDC)),
140 fESDFMD(new AliESDFMD(*esd.fESDFMD)),
141 fESDVZERO(new AliESDVZERO(*esd.fESDVZERO)),
142 fESDTZERO(new AliESDTZERO(*esd.fESDTZERO)),
143 fTPCVertex(new AliESDVertex(*esd.fTPCVertex)),
144 fSPDVertex(new AliESDVertex(*esd.fSPDVertex)),
145 fPrimaryVertex(new AliESDVertex(*esd.fPrimaryVertex)),
146 fSPDMult(new AliMultiplicity(*esd.fSPDMult)),
147 fPHOSTrigger(new AliESDCaloTrigger(*esd.fPHOSTrigger)),
148 fEMCALTrigger(new AliESDCaloTrigger(*esd.fEMCALTrigger)),
149 fESDACORDE(new AliESDACORDE(*esd.fESDACORDE)),
150 fTracks(new TClonesArray(*esd.fTracks)),
151 fMuonTracks(new TClonesArray(*esd.fMuonTracks)),
152 fPmdTracks(new TClonesArray(*esd.fPmdTracks)),
153 fTrdTracks(new TClonesArray(*esd.fTrdTracks)),
154 fV0s(new TClonesArray(*esd.fV0s)),
155 fCascades(new TClonesArray(*esd.fCascades)),
156 fKinks(new TClonesArray(*esd.fKinks)),
157 fCaloClusters(new TClonesArray(*esd.fCaloClusters)),
158 fEMCALCells(new AliESDCaloCells(*esd.fEMCALCells)),
159 fPHOSCells(new AliESDCaloCells(*esd.fPHOSCells)),
160 fErrorLogs(new TClonesArray(*esd.fErrorLogs)),
161 fESDOld(new AliESD(*esd.fESDOld)),
162 fESDFriendOld(new AliESDfriend(*esd.fESDFriendOld)),
163 fConnected(esd.fConnected),
164 fUseOwnList(esd.fUseOwnList),
165 fEMCALClusters(esd.fEMCALClusters),
166 fFirstEMCALCluster(esd.fFirstEMCALCluster),
167 fPHOSClusters(esd.fPHOSClusters),
168 fFirstPHOSCluster(esd.fFirstPHOSCluster)
171 // CKB init in the constructor list and only add here ...
176 AddObject(fESDVZERO);
177 AddObject(fESDTZERO);
178 AddObject(fTPCVertex);
179 AddObject(fSPDVertex);
180 AddObject(fPrimaryVertex);
182 AddObject(fPHOSTrigger);
183 AddObject(fEMCALTrigger);
185 AddObject(fMuonTracks);
186 AddObject(fPmdTracks);
187 AddObject(fTrdTracks);
189 AddObject(fCascades);
191 AddObject(fCaloClusters);
192 AddObject(fEMCALCells);
193 AddObject(fPHOSCells);
194 AddObject(fErrorLogs);
195 AddObject(fESDACORDE);
201 //______________________________________________________________________________
202 AliESDEvent & AliESDEvent::operator=(const AliESDEvent& source) {
204 // Assignment operator
206 if(&source == this) return *this;
207 AliVEvent::operator=(source);
209 // This assumes that the list is already created
210 // and that the virtual void Copy(Tobject&) function
211 // is correctly implemented in the derived class
212 // otherwise only TObject::Copy() will be used
215 if((fESDObjects->GetSize()==0)&&(source.fESDObjects->GetSize()>=kESDListN)){
216 // We cover the case that we do not yet have the
217 // standard content but the source has it
221 TIter next(source.GetList());
224 while ((its = next())) {
225 name.Form("%s", its->GetName());
226 TObject *mine = fESDObjects->FindObject(name.Data());
228 // not in this: can be added to list (to be implemented)
229 AliWarning(Form("%s:%d Could not find %s for copying \n",
230 (char*)__FILE__,__LINE__,name.Data()));
234 if(!its->InheritsFrom("TCollection")){
238 else if(its->InheritsFrom("TClonesArray")){
239 // Create or expand the tclonesarray pointers
240 // so we can directly copy to the object
241 TClonesArray *its_tca = (TClonesArray*)its;
242 TClonesArray *mine_tca = (TClonesArray*)mine;
244 // this leaves the capacity of the TClonesArray the same
245 // except for a factor of 2 increase when size > capacity
246 // does not release any memory occupied by the tca
247 mine_tca->ExpandCreate(its_tca->GetEntriesFast());
248 for(int i = 0;i < its_tca->GetEntriesFast();++i){
250 TObject *mine_tca_obj = mine_tca->At(i);
251 TObject *its_tca_obj = its_tca->At(i);
252 // no need to delete first
253 // pointers within the class should be handled by Copy()...
254 // Can there be Empty slots?
255 its_tca_obj->Copy(*mine_tca_obj);
259 AliWarning(Form("%s:%d cannot copy TCollection \n",
260 (char*)__FILE__,__LINE__));
264 fConnected = source.fConnected;
265 fUseOwnList = source.fUseOwnList;
266 fEMCALClusters = source.fEMCALClusters;
267 fFirstEMCALCluster = source.fFirstEMCALCluster;
268 fPHOSClusters = source.fPHOSClusters;
269 fFirstPHOSCluster = source.fFirstPHOSCluster;
277 //______________________________________________________________________________
278 AliESDEvent::~AliESDEvent()
281 // Standard destructor
284 // everthing on the list gets deleted automatically
287 if(fESDObjects&&!fConnected)
296 void AliESDEvent::Copy(TObject &obj) const {
298 // interface to TOBject::Copy
299 // Copies the content of this into obj!
300 // bascially obj = *this
302 if(this==&obj)return;
303 AliESDEvent *robj = dynamic_cast<AliESDEvent*>(&obj);
304 if(!robj)return; // not an AliESEvent
309 //______________________________________________________________________________
310 void AliESDEvent::Reset()
314 // Reset the standard contents
316 if(fESDOld)fESDOld->Reset();
317 // reset for the friends...
319 fESDFriendOld->~AliESDfriend();
320 new (fESDFriendOld) AliESDfriend();
322 // for new data we have to fetch the Pointer from the list
323 AliESDfriend *fr = (AliESDfriend*)FindListObject("AliESDfriend");
325 // delete the content
327 // make a new valid ESDfriend at the same place
328 new (fr) AliESDfriend();
331 // call reset for user supplied data?
334 void AliESDEvent::ResetStdContent()
336 // Reset the standard contents
337 if(fESDRun) fESDRun->Reset();
338 if(fHeader) fHeader->Reset();
339 if(fESDZDC) fESDZDC->Reset();
344 // reset by callin d'to /c'tor keep the pointer
345 fESDVZERO->~AliESDVZERO();
346 new (fESDVZERO) AliESDVZERO();
349 fESDACORDE->~AliESDACORDE();
350 new (fESDACORDE) AliESDACORDE();
352 if(fESDTZERO) fESDTZERO->Reset();
353 // CKB no clear/reset implemented
355 fTPCVertex->~AliESDVertex();
356 new (fTPCVertex) AliESDVertex();
357 fTPCVertex->SetName(fgkESDListName[kTPCVertex]);
360 fSPDVertex->~AliESDVertex();
361 new (fSPDVertex) AliESDVertex();
362 fSPDVertex->SetName(fgkESDListName[kSPDVertex]);
365 fPrimaryVertex->~AliESDVertex();
366 new (fPrimaryVertex) AliESDVertex();
367 fPrimaryVertex->SetName(fgkESDListName[kPrimaryVertex]);
370 fSPDMult->~AliMultiplicity();
371 new (fSPDMult) AliMultiplicity();
373 if(fPHOSTrigger)fPHOSTrigger->Reset();
374 if(fEMCALTrigger)fEMCALTrigger->Reset();
375 if(fTracks)fTracks->Delete();
376 if(fMuonTracks)fMuonTracks->Delete();
377 if(fPmdTracks)fPmdTracks->Delete();
378 if(fTrdTracks)fTrdTracks->Delete();
379 if(fV0s)fV0s->Delete();
380 if(fCascades)fCascades->Delete();
381 if(fKinks)fKinks->Delete();
382 if(fCaloClusters)fCaloClusters->Delete();
383 if(fPHOSCells)fPHOSCells->DeleteContainer();
384 if(fEMCALCells)fEMCALCells->DeleteContainer();
385 if(fErrorLogs) fErrorLogs->Delete();
387 // don't reset fconnected fConnected and the list
390 fFirstEMCALCluster=-1;
392 fFirstPHOSCluster=-1;
396 Int_t AliESDEvent::AddV0(const AliESDv0 *v) {
400 TClonesArray &fv = *fV0s;
401 Int_t idx=fV0s->GetEntriesFast();
402 new(fv[idx]) AliESDv0(*v);
406 //______________________________________________________________________________
407 void AliESDEvent::Print(Option_t *) const
410 // Print header information of the event
412 printf("ESD run information\n");
413 printf("Event # in file %d Bunch crossing # %d Orbit # %d Period # %d Run # %d Trigger %lld Magnetic field %f \n",
414 GetEventNumberInFile(),
415 GetBunchCrossNumber(),
420 GetMagneticField() );
421 printf("Vertex: (%.4f +- %.4f, %.4f +- %.4f, %.4f +- %.4f) cm\n",
422 fPrimaryVertex->GetXv(), fPrimaryVertex->GetXRes(),
423 fPrimaryVertex->GetYv(), fPrimaryVertex->GetYRes(),
424 fPrimaryVertex->GetZv(), fPrimaryVertex->GetZRes());
425 printf("Mean vertex in RUN: X=%.4f Y=%.4f cm\n",
426 GetDiamondX(),GetDiamondY());
427 printf("SPD Multiplicity. Number of tracklets %d \n",
428 fSPDMult->GetNumberOfTracklets());
429 printf("Number of tracks: \n");
430 printf(" charged %d\n", GetNumberOfTracks());
431 printf(" muon %d\n", GetNumberOfMuonTracks());
432 printf(" pmd %d\n", GetNumberOfPmdTracks());
433 printf(" trd %d\n", GetNumberOfTrdTracks());
434 printf(" v0 %d\n", GetNumberOfV0s());
435 printf(" cascades %d\n", GetNumberOfCascades());
436 printf(" kinks %d\n", GetNumberOfKinks());
437 if(fPHOSCells)printf(" PHOSCells %d\n", fPHOSCells->GetNumberOfCells());
438 else printf(" PHOSCells not in the Event\n");
439 if(fEMCALCells)printf(" EMCALCells %d\n", fEMCALCells->GetNumberOfCells());
440 else printf(" EMCALCells not in the Event\n");
441 printf(" CaloClusters %d\n", GetNumberOfCaloClusters());
442 printf(" phos %d\n", GetNumberOfPHOSClusters());
443 printf(" emcal %d\n", GetNumberOfEMCALClusters());
444 printf(" FMD %s\n", (fESDFMD ? "yes" : "no"));
445 printf(" VZERO %s\n", (fESDVZERO ? "yes" : "no"));
450 void AliESDEvent::SetESDfriend(const AliESDfriend *ev) const {
452 // Attaches the complementary info to the ESD
456 // to be sure that we set the tracks also
457 // in case of old esds
458 // if(fESDOld)CopyFromOldESD();
460 Int_t ntrk=ev->GetNumberOfTracks();
462 for (Int_t i=0; i<ntrk; i++) {
463 const AliESDfriendTrack *f=ev->GetTrack(i);
464 GetTrack(i)->SetFriendTrack(f);
468 Bool_t AliESDEvent::RemoveKink(Int_t rm) const {
469 // ---------------------------------------------------------
470 // Remove a kink candidate and references to it from ESD,
471 // if this candidate does not come from a reconstructed decay
472 // Not yet implemented...
473 // ---------------------------------------------------------
474 Int_t last=GetNumberOfKinks()-1;
475 if ((rm<0)||(rm>last)) return kFALSE;
480 Bool_t AliESDEvent::RemoveV0(Int_t rm) const {
481 // ---------------------------------------------------------
482 // Remove a V0 candidate and references to it from ESD,
483 // if this candidate does not come from a reconstructed decay
484 // ---------------------------------------------------------
485 Int_t last=GetNumberOfV0s()-1;
486 if ((rm<0)||(rm>last)) return kFALSE;
488 AliESDv0 *v0=GetV0(rm);
489 Int_t idxP=v0->GetPindex(), idxN=v0->GetNindex();
492 Int_t lastIdxP=v0->GetPindex(), lastIdxN=v0->GetNindex();
496 // Check if this V0 comes from a reconstructed decay
497 Int_t ncs=GetNumberOfCascades();
498 for (Int_t n=0; n<ncs; n++) {
499 AliESDcascade *cs=GetCascade(n);
501 Int_t csIdxP=cs->GetPindex();
502 Int_t csIdxN=cs->GetNindex();
505 if (idxN==csIdxN) return kFALSE;
507 if (csIdxP==lastIdxP)
508 if (csIdxN==lastIdxN) used++;
511 //Replace the removed V0 with the last V0
512 TClonesArray &a=*fV0s;
513 delete a.RemoveAt(rm);
515 if (rm==last) return kTRUE;
517 //v0 is pointing to the last V0 candidate...
518 new (a[rm]) AliESDv0(*v0);
519 delete a.RemoveAt(last);
521 if (!used) return kTRUE;
524 // Remap the indices of the daughters of reconstructed decays
525 for (Int_t n=0; n<ncs; n++) {
526 AliESDcascade *cs=GetCascade(n);
529 Int_t csIdxP=cs->GetPindex();
530 Int_t csIdxN=cs->GetNindex();
532 if (csIdxP==lastIdxP)
533 if (csIdxN==lastIdxN) {
534 cs->AliESDv0::SetIndex(1,idxP);
535 cs->AliESDv0::SetIndex(0,idxN);
537 if (!used) return kTRUE;
544 Bool_t AliESDEvent::RemoveTrack(Int_t rm) const {
545 // ---------------------------------------------------------
546 // Remove a track and references to it from ESD,
547 // if this track does not come from a reconstructed decay
548 // ---------------------------------------------------------
549 Int_t last=GetNumberOfTracks()-1;
550 if ((rm<0)||(rm>last)) return kFALSE;
554 // Check if this track comes from the reconstructed primary vertices
555 if (fTPCVertex && fTPCVertex->GetStatus()) {
556 UShort_t *primIdx=fTPCVertex->GetIndices();
557 Int_t n=fTPCVertex->GetNIndices();
559 Int_t idx=Int_t(primIdx[n]);
560 if (rm==idx) return kFALSE;
561 if (idx==last) used++;
564 if (fPrimaryVertex && fPrimaryVertex->GetStatus()) {
565 UShort_t *primIdx=fPrimaryVertex->GetIndices();
566 Int_t n=fPrimaryVertex->GetNIndices();
568 Int_t idx=Int_t(primIdx[n]);
569 if (rm==idx) return kFALSE;
570 if (idx==last) used++;
574 // Check if this track comes from a reconstructed decay
575 Int_t nv0=GetNumberOfV0s();
576 for (Int_t n=0; n<nv0; n++) {
577 AliESDv0 *v0=GetV0(n);
579 Int_t idx=v0->GetNindex();
580 if (rm==idx) return kFALSE;
581 if (idx==last) used++;
584 if (rm==idx) return kFALSE;
585 if (idx==last) used++;
588 Int_t ncs=GetNumberOfCascades();
589 for (Int_t n=0; n<ncs; n++) {
590 AliESDcascade *cs=GetCascade(n);
592 Int_t idx=cs->GetIndex();
593 if (rm==idx) return kFALSE;
594 if (idx==last) used++;
597 Int_t nkn=GetNumberOfKinks();
598 for (Int_t n=0; n<nkn; n++) {
599 AliESDkink *kn=GetKink(n);
601 Int_t idx=kn->GetIndex(0);
602 if (rm==idx) return kFALSE;
603 if (idx==last) used++;
606 if (rm==idx) return kFALSE;
607 if (idx==last) used++;
611 //Replace the removed track with the last track
612 TClonesArray &a=*fTracks;
613 delete a.RemoveAt(rm);
615 if (rm==last) return kTRUE;
617 AliESDtrack *t=GetTrack(last);
619 new (a[rm]) AliESDtrack(*t);
620 delete a.RemoveAt(last);
623 if (!used) return kTRUE;
626 // Remap the indices of the tracks used for the primary vertex reconstruction
627 if (fTPCVertex && fTPCVertex->GetStatus()) {
628 UShort_t *primIdx=fTPCVertex->GetIndices();
629 Int_t n=fTPCVertex->GetNIndices();
631 Int_t idx=Int_t(primIdx[n]);
633 primIdx[n]=Short_t(rm);
635 if (!used) return kTRUE;
639 if (fPrimaryVertex && fPrimaryVertex->GetStatus()) {
640 UShort_t *primIdx=fPrimaryVertex->GetIndices();
641 Int_t n=fPrimaryVertex->GetNIndices();
643 Int_t idx=Int_t(primIdx[n]);
645 primIdx[n]=Short_t(rm);
647 if (!used) return kTRUE;
652 // Remap the indices of the daughters of reconstructed decays
653 for (Int_t n=0; n<nv0; n++) {
654 AliESDv0 *v0=GetV0(n);
655 if (v0->GetIndex(0)==last) {
658 if (!used) return kTRUE;
660 if (v0->GetIndex(1)==last) {
663 if (!used) return kTRUE;
667 for (Int_t n=0; n<ncs; n++) {
668 AliESDcascade *cs=GetCascade(n);
669 if (cs->GetIndex()==last) {
672 if (!used) return kTRUE;
676 for (Int_t n=0; n<nkn; n++) {
677 AliESDkink *kn=GetKink(n);
678 if (kn->GetIndex(0)==last) {
681 if (!used) return kTRUE;
683 if (kn->GetIndex(1)==last) {
686 if (!used) return kTRUE;
694 Bool_t AliESDEvent::Clean(Float_t *cleanPars) {
696 // Remove the data which are not needed for the physics analysis.
698 // 1) Cleaning the V0 candidates
699 // ---------------------------
700 // If the cosine of the V0 pointing angle "csp" and
701 // the DCA between the daughter tracks "dca" does not satisfy
704 // csp > cleanPars[1] + dca/cleanPars[0]*(1.- cleanPars[1])
706 // an attempt to remove this V0 candidate from ESD is made.
708 // The V0 candidate gets removed if it does not belong to any
709 // recosntructed cascade decay
711 // 12.11.2007, optimal values: cleanPars[0]=0.5, cleanPars[1]=0.999
713 // 2) Cleaning the tracks
714 // ----------------------
715 // If track's transverse parameter is larger than cleanPars[2]
717 // track's longitudinal parameter is larger than cleanPars[3]
718 // an attempt to remove this track from ESD is made.
720 // The track gets removed if it does not come
721 // from a reconstructed decay
725 Float_t dcaMax=cleanPars[0];
726 Float_t cspMin=cleanPars[1];
728 Int_t nV0s=GetNumberOfV0s();
729 for (Int_t i=nV0s-1; i>=0; i--) {
730 AliESDv0 *v0=GetV0(i);
732 Float_t dca=v0->GetDcaV0Daughters();
733 Float_t csp=v0->GetV0CosineOfPointingAngle();
734 Float_t cspcut=cspMin + dca/dcaMax*(1.-cspMin);
735 if (csp > cspcut) continue;
736 if (RemoveV0(i)) rc=kTRUE;
740 Float_t dmax=cleanPars[2], zmax=cleanPars[3];
742 const AliESDVertex *vertex=GetPrimaryVertexSPD();
743 Bool_t vtxOK=vertex->GetStatus();
745 Int_t nTracks=GetNumberOfTracks();
746 for (Int_t i=nTracks-1; i>=0; i--) {
747 AliESDtrack *track=GetTrack(i);
748 Float_t xy,z; track->GetImpactParameters(xy,z);
749 if ((TMath::Abs(xy) > dmax) || (vtxOK && (TMath::Abs(z) > zmax))) {
750 if (RemoveTrack(i)) rc=kTRUE;
757 Int_t AliESDEvent::AddTrack(const AliESDtrack *t)
760 TClonesArray &ftr = *fTracks;
761 AliESDtrack * track = new(ftr[fTracks->GetEntriesFast()])AliESDtrack(*t);
762 track->SetID(fTracks->GetEntriesFast()-1);
763 return track->GetID();
766 void AliESDEvent::AddMuonTrack(const AliESDMuonTrack *t)
768 TClonesArray &fmu = *fMuonTracks;
769 new(fmu[fMuonTracks->GetEntriesFast()]) AliESDMuonTrack(*t);
772 void AliESDEvent::AddPmdTrack(const AliESDPmdTrack *t)
774 TClonesArray &fpmd = *fPmdTracks;
775 new(fpmd[fPmdTracks->GetEntriesFast()]) AliESDPmdTrack(*t);
778 void AliESDEvent::AddTrdTrack(const AliESDTrdTrack *t)
780 TClonesArray &ftrd = *fTrdTracks;
781 new(ftrd[fTrdTracks->GetEntriesFast()]) AliESDTrdTrack(*t);
787 Int_t AliESDEvent::AddKink(const AliESDkink *c)
790 TClonesArray &fk = *fKinks;
791 AliESDkink * kink = new(fk[fKinks->GetEntriesFast()]) AliESDkink(*c);
792 kink->SetID(fKinks->GetEntriesFast()); // CKB different from the other imps..
793 return fKinks->GetEntriesFast()-1;
797 void AliESDEvent::AddCascade(const AliESDcascade *c)
799 TClonesArray &fc = *fCascades;
800 new(fc[fCascades->GetEntriesFast()]) AliESDcascade(*c);
804 Int_t AliESDEvent::AddCaloCluster(const AliESDCaloCluster *c)
807 TClonesArray &fc = *fCaloClusters;
808 AliESDCaloCluster *clus = new(fc[fCaloClusters->GetEntriesFast()]) AliESDCaloCluster(*c);
809 clus->SetID(fCaloClusters->GetEntriesFast()-1);
810 return fCaloClusters->GetEntriesFast()-1;
814 void AliESDEvent::AddRawDataErrorLog(const AliRawDataErrorLog *log) const {
815 TClonesArray &errlogs = *fErrorLogs;
816 new(errlogs[errlogs.GetEntriesFast()]) AliRawDataErrorLog(*log);
819 void AliESDEvent::SetPrimaryVertexTPC(const AliESDVertex *vertex)
821 // Set the TPC vertex
822 // use already allocated space
824 *fTPCVertex = *vertex;
825 fTPCVertex->SetName(fgkESDListName[kTPCVertex]);
829 void AliESDEvent::SetPrimaryVertexSPD(const AliESDVertex *vertex)
831 // Set the SPD vertex
832 // use already allocated space
834 *fSPDVertex = *vertex;
835 fSPDVertex->SetName(fgkESDListName[kSPDVertex]);
839 void AliESDEvent::SetPrimaryVertex(const AliESDVertex *vertex)
841 // Set the primary vertex
842 // use already allocated space
844 *fPrimaryVertex = *vertex;
845 fPrimaryVertex->SetName(fgkESDListName[kPrimaryVertex]);
849 void AliESDEvent::SetMultiplicity(const AliMultiplicity *mul)
851 // Set the SPD Multiplicity
858 void AliESDEvent::SetFMDData(AliESDFMD * obj)
860 // use already allocated space
866 void AliESDEvent::SetVZEROData(AliESDVZERO * obj)
868 // use already allocated space
873 void AliESDEvent::SetACORDEData(AliESDACORDE * obj)
880 void AliESDEvent::GetESDfriend(AliESDfriend *ev) const
883 // Extracts the complementary info from the ESD
887 Int_t ntrk=GetNumberOfTracks();
889 for (Int_t i=0; i<ntrk; i++) {
890 AliESDtrack *t=GetTrack(i);
891 const AliESDfriendTrack *f=t->GetFriendTrack();
894 t->ReleaseESDfriendTrack();// Not to have two copies of "friendTrack"
898 AliESDfriend *fr = (AliESDfriend*)(const_cast<AliESDEvent*>(this)->FindListObject("AliESDfriend"));
899 if (fr) ev->SetVZEROfriend(fr->GetVZEROfriend());
902 void AliESDEvent::AddObject(TObject* obj)
904 // Add an object to the list of object.
905 // Please be aware that in order to increase performance you should
906 // refrain from using TObjArrays (if possible). Use TClonesArrays, instead.
907 fESDObjects->SetOwner(kTRUE);
908 fESDObjects->AddLast(obj);
912 void AliESDEvent::GetStdContent()
914 // set pointers for standard content
915 // get by name much safer and not a big overhead since not called very often
917 fESDRun = (AliESDRun*)fESDObjects->FindObject(fgkESDListName[kESDRun]);
918 fHeader = (AliESDHeader*)fESDObjects->FindObject(fgkESDListName[kHeader]);
919 fESDZDC = (AliESDZDC*)fESDObjects->FindObject(fgkESDListName[kESDZDC]);
920 fESDFMD = (AliESDFMD*)fESDObjects->FindObject(fgkESDListName[kESDFMD]);
921 fESDVZERO = (AliESDVZERO*)fESDObjects->FindObject(fgkESDListName[kESDVZERO]);
922 fESDTZERO = (AliESDTZERO*)fESDObjects->FindObject(fgkESDListName[kESDTZERO]);
923 fTPCVertex = (AliESDVertex*)fESDObjects->FindObject(fgkESDListName[kTPCVertex]);
924 fSPDVertex = (AliESDVertex*)fESDObjects->FindObject(fgkESDListName[kSPDVertex]);
925 fPrimaryVertex = (AliESDVertex*)fESDObjects->FindObject(fgkESDListName[kPrimaryVertex]);
926 fSPDMult = (AliMultiplicity*)fESDObjects->FindObject(fgkESDListName[kSPDMult]);
927 fPHOSTrigger = (AliESDCaloTrigger*)fESDObjects->FindObject(fgkESDListName[kPHOSTrigger]);
928 fEMCALTrigger = (AliESDCaloTrigger*)fESDObjects->FindObject(fgkESDListName[kEMCALTrigger]);
929 fTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kTracks]);
930 fMuonTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kMuonTracks]);
931 fPmdTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kPmdTracks]);
932 fTrdTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kTrdTracks]);
933 fV0s = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kV0s]);
934 fCascades = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kCascades]);
935 fKinks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kKinks]);
936 fCaloClusters = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kCaloClusters]);
937 fEMCALCells = (AliESDCaloCells*)fESDObjects->FindObject(fgkESDListName[kEMCALCells]);
938 fPHOSCells = (AliESDCaloCells*)fESDObjects->FindObject(fgkESDListName[kPHOSCells]);
939 fErrorLogs = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kErrorLogs]);
940 fESDACORDE = (AliESDACORDE*)fESDObjects->FindObject(fgkESDListName[kESDACORDE]);
944 void AliESDEvent::SetStdNames(){
945 // Set the names of the standard contents
947 if(fESDObjects->GetEntries()==kESDListN){
948 for(int i = 0;i < fESDObjects->GetEntries();i++){
949 TObject *fObj = fESDObjects->At(i);
950 if(fObj->InheritsFrom("TNamed")){
951 ((TNamed*)fObj)->SetName(fgkESDListName[i]);
953 else if(fObj->InheritsFrom("TClonesArray")){
954 ((TClonesArray*)fObj)->SetName(fgkESDListName[i]);
959 printf("%s:%d SetStdNames() Wrong number of Std Entries \n",(char*)__FILE__,__LINE__);
964 void AliESDEvent::CreateStdContent(Bool_t bUseThisList){
965 fUseOwnList = bUseThisList;
969 void AliESDEvent::CreateStdContent()
971 // create the standard AOD content and set pointers
973 // create standard objects and add them to the TList of objects
974 AddObject(new AliESDRun());
975 AddObject(new AliESDHeader());
976 AddObject(new AliESDZDC());
977 AddObject(new AliESDFMD());
978 AddObject(new AliESDVZERO());
979 AddObject(new AliESDTZERO());
980 AddObject(new AliESDVertex());
981 AddObject(new AliESDVertex());
982 AddObject(new AliESDVertex());
983 AddObject(new AliMultiplicity());
984 AddObject(new AliESDCaloTrigger());
985 AddObject(new AliESDCaloTrigger());
986 AddObject(new TClonesArray("AliESDtrack",0));
987 AddObject(new TClonesArray("AliESDMuonTrack",0));
988 AddObject(new TClonesArray("AliESDPmdTrack",0));
989 AddObject(new TClonesArray("AliESDTrdTrack",0));
990 AddObject(new TClonesArray("AliESDv0",0));
991 AddObject(new TClonesArray("AliESDcascade",0));
992 AddObject(new TClonesArray("AliESDkink",0));
993 AddObject(new TClonesArray("AliESDCaloCluster",0));
994 AddObject(new AliESDCaloCells());
995 AddObject(new AliESDCaloCells());
996 AddObject(new TClonesArray("AliRawDataErrorLog",0));
997 AddObject(new AliESDACORDE());
999 // check the order of the indices against enum...
1003 // read back pointers
1007 TObject* AliESDEvent::FindListObject(const char *name){
1009 // Find object with name "name" in the list of branches
1012 return fESDObjects->FindObject(name);
1017 Int_t AliESDEvent::GetPHOSClusters(TRefArray *clusters) const
1019 // fills the provided TRefArray with all found phos clusters
1023 AliESDCaloCluster *cl = 0;
1024 for (Int_t i = 0; i < GetNumberOfCaloClusters(); i++) {
1026 if ( (cl = GetCaloCluster(i)) ) {
1029 AliDebug(1,Form("IsPHOS cluster %d Size: %d \n",i,clusters->GetEntriesFast()));
1033 return clusters->GetEntriesFast();
1036 Int_t AliESDEvent::GetEMCALClusters(TRefArray *clusters) const
1038 // fills the provided TRefArray with all found emcal clusters
1042 AliESDCaloCluster *cl = 0;
1043 for (Int_t i = 0; i < GetNumberOfCaloClusters(); i++) {
1045 if ( (cl = GetCaloCluster(i)) ) {
1048 AliDebug(1,Form("IsEMCAL cluster %d Size: %d \n",i,clusters->GetEntriesFast()));
1052 return clusters->GetEntriesFast();
1055 const void AliESDEvent::WriteToTree(TTree* tree) const {
1056 // Book the branches as in TTree::Branch(TCollection*)
1057 // but add a "." at the end of top level branches which are
1058 // not a TClonesArray
1062 TIter next(fESDObjects);
1063 const Int_t kSplitlevel = 99; // default value in TTree::Branch()
1064 const Int_t kBufsize = 32000; // default value in TTree::Branch()
1067 while ((obj = next())) {
1068 branchname.Form("%s", obj->GetName());
1069 if ((kSplitlevel > 1) && !obj->InheritsFrom(TClonesArray::Class())) {
1070 if(!branchname.EndsWith("."))branchname += ".";
1072 tree->Bronch(branchname, obj->ClassName(), fESDObjects->GetObjectRef(obj),
1073 kBufsize, kSplitlevel - 1);
1079 void AliESDEvent::ReadFromTree(TTree *tree, Option_t* /*opt*/){
1081 // Connect the ESDEvent to a tree
1084 Printf("%s %d AliESDEvent::ReadFromTree() Zero Pointer to Tree \n",(char*)__FILE__,__LINE__);
1088 if(!tree->GetTree())tree->LoadTree(0);
1090 // if we find the "ESD" branch on the tree we do have the old structure
1091 if(tree->GetBranch("ESD")) {
1092 char ** address = (char **)(tree->GetBranch("ESD")->GetAddress());
1093 // do we have the friend branch
1094 TBranch * esdFB = tree->GetBranch("ESDfriend.");
1095 char ** addressF = 0;
1096 if(esdFB)addressF = (char **)(esdFB->GetAddress());
1098 printf("%s %d AliESDEvent::ReadFromTree() Reading old Tree \n",(char*)__FILE__,__LINE__);
1099 tree->SetBranchAddress("ESD", &fESDOld);
1101 tree->SetBranchAddress("ESDfriend.",&fESDFriendOld);
1104 printf("%s %d AliESDEvent::ReadFromTree() Reading old Tree \n",(char*)__FILE__,__LINE__);
1105 printf("%s %d Branch already connected. Using existing branch address. \n",(char*)__FILE__,__LINE__);
1106 fESDOld = (AliESD*) (*address);
1107 // addressF can still be 0, since branch needs to switched on
1108 if(addressF)fESDFriendOld = (AliESDfriend*) (*addressF);
1111 // have already connected the old ESD structure... ?
1112 // reuse also the pointer of the AlliESDEvent
1113 // otherwise create new ones
1114 TList* connectedList = (TList*) (tree->GetUserInfo()->FindObject("ESDObjectsConnectedToTree"));
1117 // If connected use the connected list of objects
1118 if(fESDObjects!= connectedList){
1119 // protect when called twice
1120 fESDObjects->Delete();
1121 fESDObjects = connectedList;
1126 // The pointer to the friend changes when called twice via InitIO
1127 // since AliESDEvent is deleted
1128 TObject* oldf = FindListObject("AliESDfriend");
1131 newf = (TObject*)*addressF;
1133 if(newf!=0&&oldf!=newf){
1134 // remove the old reference
1135 // Should we also delete it? Or is this handled in TTree I/O
1136 // since it is created by the first SetBranchAddress
1137 fESDObjects->Remove(oldf);
1139 fESDObjects->Add(newf);
1146 CreateStdContent(); // create for copy
1147 // if we have the esdfriend add it, so we always can access it via the userinfo
1148 if(fESDFriendOld)AddObject(fESDFriendOld);
1149 // we are not owner of the list objects
1150 // must not delete it
1151 fESDObjects->SetOwner(kFALSE);
1152 fESDObjects->SetName("ESDObjectsConnectedToTree");
1153 tree->GetUserInfo()->Add(fESDObjects);
1161 // Try to find AliESDEvent
1162 AliESDEvent *esdEvent = 0;
1163 esdEvent = (AliESDEvent*)tree->GetTree()->GetUserInfo()->FindObject("AliESDEvent");
1165 // Check if already connected to tree
1166 TList* connectedList = (TList*) (tree->GetUserInfo()->FindObject("ESDObjectsConnectedToTree"));
1167 if (connectedList) {
1168 // If connected use the connected list if objects
1169 fESDObjects->Delete();
1170 fESDObjects = connectedList;
1177 // prevent a memory leak when reading back the TList
1182 // create a new TList from the UserInfo TList...
1183 // copy constructor does not work...
1184 fESDObjects = (TList*)(esdEvent->GetList()->Clone());
1185 fESDObjects->SetOwner(kFALSE);
1187 else if ( fESDObjects->GetEntries()==0){
1188 // at least create the std content if we want to read to our list
1193 // we only need new things in the list if we do no already have it..
1194 // TODO just add new entries
1196 if(fESDObjects->GetEntries()<kESDListN){
1197 printf("%s %d AliESDEvent::ReadFromTree() TList contains less than the standard contents %d < %d \n",
1198 (char*)__FILE__,__LINE__,fESDObjects->GetEntries(),kESDListN);
1200 // set the branch addresses
1201 TIter next(fESDObjects);
1203 while((el=(TNamed*)next())){
1204 TString bname(el->GetName());
1205 if(bname.CompareTo("AliESDfriend")==0)
1207 // AliESDfriend does not have a name ...
1208 tree->SetBranchAddress("ESDfriend.",fESDObjects->GetObjectRef(el));
1211 // check if branch exists under this Name
1212 TBranch *br = tree->GetBranch(bname.Data());
1214 tree->SetBranchAddress(bname.Data(),fESDObjects->GetObjectRef(el));
1217 br = tree->GetBranch(Form("%s.",bname.Data()));
1219 tree->SetBranchAddress(Form("%s.",bname.Data()),fESDObjects->GetObjectRef(el));
1222 printf("%s %d AliESDEvent::ReadFromTree() No Branch found with Name %s or %s. \n",
1223 (char*)__FILE__,__LINE__,bname.Data(),bname.Data());
1230 // when reading back we are not owner of the list
1231 // must not delete it
1232 fESDObjects->SetOwner(kFALSE);
1233 fESDObjects->SetName("ESDObjectsConnectedToTree");
1234 // we are not owner of the list objects
1235 // must not delete it
1236 tree->GetUserInfo()->Add(fESDObjects);
1240 // we can't get the list from the user data, create standard content
1241 // and set it by hand (no ESDfriend at the moment
1243 TIter next(fESDObjects);
1245 while((el=(TNamed*)next())){
1246 TString bname(el->GetName());
1247 TBranch *br = tree->GetBranch(bname.Data());
1249 tree->SetBranchAddress(bname.Data(),fESDObjects->GetObjectRef(el));
1252 br = tree->GetBranch(Form("%s.",bname.Data()));
1254 tree->SetBranchAddress(Form("%s.",bname.Data()),fESDObjects->GetObjectRef(el));
1259 // when reading back we are not owner of the list
1260 // must not delete it
1261 fESDObjects->SetOwner(kFALSE);
1266 void AliESDEvent::CopyFromOldESD()
1268 // Method which copies over everthing from the old esd structure to the
1273 SetRunNumber(fESDOld->GetRunNumber());
1274 SetPeriodNumber(fESDOld->GetPeriodNumber());
1275 SetMagneticField(fESDOld->GetMagneticField());
1277 // leave out diamond ...
1278 // SetDiamond(const AliESDVertex *vertex) { fESDRun->SetDiamond(vertex);}
1281 SetTriggerMask(fESDOld->GetTriggerMask());
1282 SetOrbitNumber(fESDOld->GetOrbitNumber());
1283 SetTimeStamp(fESDOld->GetTimeStamp());
1284 SetEventType(fESDOld->GetEventType());
1285 SetEventNumberInFile(fESDOld->GetEventNumberInFile());
1286 SetBunchCrossNumber(fESDOld->GetBunchCrossNumber());
1287 SetTriggerCluster(fESDOld->GetTriggerCluster());
1291 SetZDC(fESDOld->GetZDCN1Energy(),
1292 fESDOld->GetZDCP1Energy(),
1293 fESDOld->GetZDCEMEnergy(),
1295 fESDOld->GetZDCN2Energy(),
1296 fESDOld->GetZDCP2Energy(),
1297 fESDOld->GetZDCParticipants());
1301 if(fESDOld->GetFMDData())SetFMDData(fESDOld->GetFMDData());
1305 SetT0zVertex(fESDOld->GetT0zVertex());
1306 SetT0(fESDOld->GetT0());
1310 if (fESDOld->GetVZEROData()) SetVZEROData(fESDOld->GetVZEROData());
1312 if(fESDOld->GetVertex())SetPrimaryVertexSPD(fESDOld->GetVertex());
1314 if(fESDOld->GetPrimaryVertex())SetPrimaryVertex(fESDOld->GetPrimaryVertex());
1316 if(fESDOld->GetMultiplicity())SetMultiplicity(fESDOld->GetMultiplicity());
1318 for(int i = 0;i<fESDOld->GetNumberOfTracks();i++){
1319 AddTrack(fESDOld->GetTrack(i));
1322 for(int i = 0;i<fESDOld->GetNumberOfMuonTracks();i++){
1323 AddMuonTrack(fESDOld->GetMuonTrack(i));
1326 for(int i = 0;i<fESDOld->GetNumberOfPmdTracks();i++){
1327 AddPmdTrack(fESDOld->GetPmdTrack(i));
1330 for(int i = 0;i<fESDOld->GetNumberOfTrdTracks();i++){
1331 AddTrdTrack(fESDOld->GetTrdTrack(i));
1334 for(int i = 0;i<fESDOld->GetNumberOfV0s();i++){
1335 AddV0(fESDOld->GetV0(i));
1338 for(int i = 0;i<fESDOld->GetNumberOfCascades();i++){
1339 AddCascade(fESDOld->GetCascade(i));
1342 for(int i = 0;i<fESDOld->GetNumberOfKinks();i++){
1343 AddKink(fESDOld->GetKink(i));
1347 for(int i = 0;i<fESDOld->GetNumberOfCaloClusters();i++){
1348 AddCaloCluster(fESDOld->GetCaloCluster(i));