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
216 if((fESDObjects->GetSize()==0)&&(source.fESDObjects->GetSize()>=kESDListN)){
217 // We cover the case that we do not yet have the
218 // standard content but the source has it
222 TIter next(source.GetList());
225 while ((its = next())) {
226 name.Form("%s", its->GetName());
227 TObject *mine = fESDObjects->FindObject(name.Data());
229 TClass* pClass=TClass::GetClass(its->ClassName());
231 AliWarning(Form("Can not find class description for entry %s (%s)\n",
232 its->ClassName(), name.Data()));
236 mine=(TObject*)pClass->New();
238 // not in this: can be added to list
239 AliWarning(Form("%s:%d Could not find %s for copying \n",
240 (char*)__FILE__,__LINE__,name.Data()));
243 if(mine->InheritsFrom("TNamed")){
244 ((TNamed*)mine)->SetName(name);
246 else if(mine->InheritsFrom("TCollection")){
247 if(mine->InheritsFrom("TClonesArray"))
248 dynamic_cast<TClonesArray*>(mine)->SetClass(dynamic_cast<TClonesArray*>(its)->GetClass());
249 dynamic_cast<TCollection*>(mine)->SetName(name);
251 AliDebug(1, Form("adding object %s of type %s", mine->GetName(), mine->ClassName()));
255 if(!its->InheritsFrom("TCollection")){
259 else if(its->InheritsFrom("TClonesArray")){
260 // Create or expand the tclonesarray pointers
261 // so we can directly copy to the object
262 TClonesArray *its_tca = (TClonesArray*)its;
263 TClonesArray *mine_tca = (TClonesArray*)mine;
265 // this leaves the capacity of the TClonesArray the same
266 // except for a factor of 2 increase when size > capacity
267 // does not release any memory occupied by the tca
268 mine_tca->ExpandCreate(its_tca->GetEntriesFast());
269 for(int i = 0;i < its_tca->GetEntriesFast();++i){
271 TObject *mine_tca_obj = mine_tca->At(i);
272 TObject *its_tca_obj = its_tca->At(i);
273 // no need to delete first
274 // pointers within the class should be handled by Copy()...
275 // Can there be Empty slots?
276 its_tca_obj->Copy(*mine_tca_obj);
280 AliWarning(Form("%s:%d cannot copy TCollection \n",
281 (char*)__FILE__,__LINE__));
285 fConnected = source.fConnected;
286 fUseOwnList = source.fUseOwnList;
287 fEMCALClusters = source.fEMCALClusters;
288 fFirstEMCALCluster = source.fFirstEMCALCluster;
289 fPHOSClusters = source.fPHOSClusters;
290 fFirstPHOSCluster = source.fFirstPHOSCluster;
298 //______________________________________________________________________________
299 AliESDEvent::~AliESDEvent()
302 // Standard destructor
305 // everthing on the list gets deleted automatically
308 if(fESDObjects&&!fConnected)
317 void AliESDEvent::Copy(TObject &obj) const {
319 // interface to TOBject::Copy
320 // Copies the content of this into obj!
321 // bascially obj = *this
323 if(this==&obj)return;
324 AliESDEvent *robj = dynamic_cast<AliESDEvent*>(&obj);
325 if(!robj)return; // not an AliESEvent
330 //______________________________________________________________________________
331 void AliESDEvent::Reset()
335 // Std content + Non std content
337 // Reset the standard contents
340 // reset for the old data without AliESDEvent...
341 if(fESDOld)fESDOld->Reset();
343 fESDFriendOld->~AliESDfriend();
344 new (fESDFriendOld) AliESDfriend();
348 if(fESDObjects->GetSize()>kESDListN){
349 // we have non std content
350 // this also covers esdfriends
351 for(int i = kESDListN;i < fESDObjects->GetSize();++i){
352 TObject *pObject = fESDObjects->At(i);
354 if(pObject->InheritsFrom(TClonesArray::Class())){
355 ((TClonesArray*)pObject)->Delete();
357 else if(!pObject->InheritsFrom(TCollection::Class())){
358 ResetWithPlacementNew(pObject);
361 AliWarning(Form("No reset for %s (%s)\n",
362 pObject->ClassName()));
369 Bool_t AliESDEvent::ResetWithPlacementNew(TObject *pObject){
370 Long_t dtoronly = TObject::GetDtorOnly();
371 TClass *pClass = TClass::GetClass(pObject->ClassName());
372 TObject::SetDtorOnly(pObject);
374 // Recreate with placement new
375 pClass->New(pObject);
376 // Restore the state.
377 TObject::SetDtorOnly((void*)dtoronly);
381 void AliESDEvent::ResetStdContent()
383 // Reset the standard contents
384 if(fESDRun) fESDRun->Reset();
385 if(fHeader) fHeader->Reset();
386 if(fESDZDC) fESDZDC->Reset();
391 // reset by callin d'to /c'tor keep the pointer
392 fESDVZERO->~AliESDVZERO();
393 new (fESDVZERO) AliESDVZERO();
396 fESDACORDE->~AliESDACORDE();
397 new (fESDACORDE) AliESDACORDE();
399 if(fESDTZERO) fESDTZERO->Reset();
400 // CKB no clear/reset implemented
402 fTPCVertex->~AliESDVertex();
403 new (fTPCVertex) AliESDVertex();
404 fTPCVertex->SetName(fgkESDListName[kTPCVertex]);
407 fSPDVertex->~AliESDVertex();
408 new (fSPDVertex) AliESDVertex();
409 fSPDVertex->SetName(fgkESDListName[kSPDVertex]);
412 fPrimaryVertex->~AliESDVertex();
413 new (fPrimaryVertex) AliESDVertex();
414 fPrimaryVertex->SetName(fgkESDListName[kPrimaryVertex]);
417 fSPDMult->~AliMultiplicity();
418 new (fSPDMult) AliMultiplicity();
420 if(fPHOSTrigger)fPHOSTrigger->Reset();
421 if(fEMCALTrigger)fEMCALTrigger->Reset();
422 if(fTracks)fTracks->Delete();
423 if(fMuonTracks)fMuonTracks->Delete();
424 if(fPmdTracks)fPmdTracks->Delete();
425 if(fTrdTracks)fTrdTracks->Delete();
426 if(fV0s)fV0s->Delete();
427 if(fCascades)fCascades->Delete();
428 if(fKinks)fKinks->Delete();
429 if(fCaloClusters)fCaloClusters->Delete();
430 if(fPHOSCells)fPHOSCells->DeleteContainer();
431 if(fEMCALCells)fEMCALCells->DeleteContainer();
432 if(fErrorLogs) fErrorLogs->Delete();
434 // don't reset fconnected fConnected and the list
437 fFirstEMCALCluster=-1;
439 fFirstPHOSCluster=-1;
443 Int_t AliESDEvent::AddV0(const AliESDv0 *v) {
447 TClonesArray &fv = *fV0s;
448 Int_t idx=fV0s->GetEntriesFast();
449 new(fv[idx]) AliESDv0(*v);
453 //______________________________________________________________________________
454 void AliESDEvent::Print(Option_t *) const
457 // Print header information of the event
459 printf("ESD run information\n");
460 printf("Event # in file %d Bunch crossing # %d Orbit # %d Period # %d Run # %d Trigger %lld Magnetic field %f \n",
461 GetEventNumberInFile(),
462 GetBunchCrossNumber(),
467 GetMagneticField() );
468 printf("Vertex: (%.4f +- %.4f, %.4f +- %.4f, %.4f +- %.4f) cm\n",
469 fPrimaryVertex->GetXv(), fPrimaryVertex->GetXRes(),
470 fPrimaryVertex->GetYv(), fPrimaryVertex->GetYRes(),
471 fPrimaryVertex->GetZv(), fPrimaryVertex->GetZRes());
472 printf("Mean vertex in RUN: X=%.4f Y=%.4f cm\n",
473 GetDiamondX(),GetDiamondY());
474 printf("SPD Multiplicity. Number of tracklets %d \n",
475 fSPDMult->GetNumberOfTracklets());
476 printf("Number of tracks: \n");
477 printf(" charged %d\n", GetNumberOfTracks());
478 printf(" muon %d\n", GetNumberOfMuonTracks());
479 printf(" pmd %d\n", GetNumberOfPmdTracks());
480 printf(" trd %d\n", GetNumberOfTrdTracks());
481 printf(" v0 %d\n", GetNumberOfV0s());
482 printf(" cascades %d\n", GetNumberOfCascades());
483 printf(" kinks %d\n", GetNumberOfKinks());
484 if(fPHOSCells)printf(" PHOSCells %d\n", fPHOSCells->GetNumberOfCells());
485 else printf(" PHOSCells not in the Event\n");
486 if(fEMCALCells)printf(" EMCALCells %d\n", fEMCALCells->GetNumberOfCells());
487 else printf(" EMCALCells not in the Event\n");
488 printf(" CaloClusters %d\n", GetNumberOfCaloClusters());
489 printf(" phos %d\n", GetNumberOfPHOSClusters());
490 printf(" emcal %d\n", GetNumberOfEMCALClusters());
491 printf(" FMD %s\n", (fESDFMD ? "yes" : "no"));
492 printf(" VZERO %s\n", (fESDVZERO ? "yes" : "no"));
497 void AliESDEvent::SetESDfriend(const AliESDfriend *ev) const {
499 // Attaches the complementary info to the ESD
503 // to be sure that we set the tracks also
504 // in case of old esds
505 // if(fESDOld)CopyFromOldESD();
507 Int_t ntrk=ev->GetNumberOfTracks();
509 for (Int_t i=0; i<ntrk; i++) {
510 const AliESDfriendTrack *f=ev->GetTrack(i);
511 GetTrack(i)->SetFriendTrack(f);
515 Bool_t AliESDEvent::RemoveKink(Int_t rm) const {
516 // ---------------------------------------------------------
517 // Remove a kink candidate and references to it from ESD,
518 // if this candidate does not come from a reconstructed decay
519 // Not yet implemented...
520 // ---------------------------------------------------------
521 Int_t last=GetNumberOfKinks()-1;
522 if ((rm<0)||(rm>last)) return kFALSE;
527 Bool_t AliESDEvent::RemoveV0(Int_t rm) const {
528 // ---------------------------------------------------------
529 // Remove a V0 candidate and references to it from ESD,
530 // if this candidate does not come from a reconstructed decay
531 // ---------------------------------------------------------
532 Int_t last=GetNumberOfV0s()-1;
533 if ((rm<0)||(rm>last)) return kFALSE;
535 AliESDv0 *v0=GetV0(rm);
536 Int_t idxP=v0->GetPindex(), idxN=v0->GetNindex();
539 Int_t lastIdxP=v0->GetPindex(), lastIdxN=v0->GetNindex();
543 // Check if this V0 comes from a reconstructed decay
544 Int_t ncs=GetNumberOfCascades();
545 for (Int_t n=0; n<ncs; n++) {
546 AliESDcascade *cs=GetCascade(n);
548 Int_t csIdxP=cs->GetPindex();
549 Int_t csIdxN=cs->GetNindex();
552 if (idxN==csIdxN) return kFALSE;
554 if (csIdxP==lastIdxP)
555 if (csIdxN==lastIdxN) used++;
558 //Replace the removed V0 with the last V0
559 TClonesArray &a=*fV0s;
560 delete a.RemoveAt(rm);
562 if (rm==last) return kTRUE;
564 //v0 is pointing to the last V0 candidate...
565 new (a[rm]) AliESDv0(*v0);
566 delete a.RemoveAt(last);
568 if (!used) return kTRUE;
571 // Remap the indices of the daughters of reconstructed decays
572 for (Int_t n=0; n<ncs; n++) {
573 AliESDcascade *cs=GetCascade(n);
576 Int_t csIdxP=cs->GetPindex();
577 Int_t csIdxN=cs->GetNindex();
579 if (csIdxP==lastIdxP)
580 if (csIdxN==lastIdxN) {
581 cs->AliESDv0::SetIndex(1,idxP);
582 cs->AliESDv0::SetIndex(0,idxN);
584 if (!used) return kTRUE;
591 Bool_t AliESDEvent::RemoveTrack(Int_t rm) const {
592 // ---------------------------------------------------------
593 // Remove a track and references to it from ESD,
594 // if this track does not come from a reconstructed decay
595 // ---------------------------------------------------------
596 Int_t last=GetNumberOfTracks()-1;
597 if ((rm<0)||(rm>last)) return kFALSE;
601 // Check if this track comes from the reconstructed primary vertices
602 if (fTPCVertex && fTPCVertex->GetStatus()) {
603 UShort_t *primIdx=fTPCVertex->GetIndices();
604 Int_t n=fTPCVertex->GetNIndices();
606 Int_t idx=Int_t(primIdx[n]);
607 if (rm==idx) return kFALSE;
608 if (idx==last) used++;
611 if (fPrimaryVertex && fPrimaryVertex->GetStatus()) {
612 UShort_t *primIdx=fPrimaryVertex->GetIndices();
613 Int_t n=fPrimaryVertex->GetNIndices();
615 Int_t idx=Int_t(primIdx[n]);
616 if (rm==idx) return kFALSE;
617 if (idx==last) used++;
621 // Check if this track comes from a reconstructed decay
622 Int_t nv0=GetNumberOfV0s();
623 for (Int_t n=0; n<nv0; n++) {
624 AliESDv0 *v0=GetV0(n);
626 Int_t idx=v0->GetNindex();
627 if (rm==idx) return kFALSE;
628 if (idx==last) used++;
631 if (rm==idx) return kFALSE;
632 if (idx==last) used++;
635 Int_t ncs=GetNumberOfCascades();
636 for (Int_t n=0; n<ncs; n++) {
637 AliESDcascade *cs=GetCascade(n);
639 Int_t idx=cs->GetIndex();
640 if (rm==idx) return kFALSE;
641 if (idx==last) used++;
644 Int_t nkn=GetNumberOfKinks();
645 for (Int_t n=0; n<nkn; n++) {
646 AliESDkink *kn=GetKink(n);
648 Int_t idx=kn->GetIndex(0);
649 if (rm==idx) return kFALSE;
650 if (idx==last) used++;
653 if (rm==idx) return kFALSE;
654 if (idx==last) used++;
658 //Replace the removed track with the last track
659 TClonesArray &a=*fTracks;
660 delete a.RemoveAt(rm);
662 if (rm==last) return kTRUE;
664 AliESDtrack *t=GetTrack(last);
666 new (a[rm]) AliESDtrack(*t);
667 delete a.RemoveAt(last);
670 if (!used) return kTRUE;
673 // Remap the indices of the tracks used for the primary vertex reconstruction
674 if (fTPCVertex && fTPCVertex->GetStatus()) {
675 UShort_t *primIdx=fTPCVertex->GetIndices();
676 Int_t n=fTPCVertex->GetNIndices();
678 Int_t idx=Int_t(primIdx[n]);
680 primIdx[n]=Short_t(rm);
682 if (!used) return kTRUE;
686 if (fPrimaryVertex && fPrimaryVertex->GetStatus()) {
687 UShort_t *primIdx=fPrimaryVertex->GetIndices();
688 Int_t n=fPrimaryVertex->GetNIndices();
690 Int_t idx=Int_t(primIdx[n]);
692 primIdx[n]=Short_t(rm);
694 if (!used) return kTRUE;
699 // Remap the indices of the daughters of reconstructed decays
700 for (Int_t n=0; n<nv0; n++) {
701 AliESDv0 *v0=GetV0(n);
702 if (v0->GetIndex(0)==last) {
705 if (!used) return kTRUE;
707 if (v0->GetIndex(1)==last) {
710 if (!used) return kTRUE;
714 for (Int_t n=0; n<ncs; n++) {
715 AliESDcascade *cs=GetCascade(n);
716 if (cs->GetIndex()==last) {
719 if (!used) return kTRUE;
723 for (Int_t n=0; n<nkn; n++) {
724 AliESDkink *kn=GetKink(n);
725 if (kn->GetIndex(0)==last) {
728 if (!used) return kTRUE;
730 if (kn->GetIndex(1)==last) {
733 if (!used) return kTRUE;
741 Bool_t AliESDEvent::Clean(Float_t *cleanPars) {
743 // Remove the data which are not needed for the physics analysis.
745 // 1) Cleaning the V0 candidates
746 // ---------------------------
747 // If the cosine of the V0 pointing angle "csp" and
748 // the DCA between the daughter tracks "dca" does not satisfy
751 // csp > cleanPars[1] + dca/cleanPars[0]*(1.- cleanPars[1])
753 // an attempt to remove this V0 candidate from ESD is made.
755 // The V0 candidate gets removed if it does not belong to any
756 // recosntructed cascade decay
758 // 12.11.2007, optimal values: cleanPars[0]=0.5, cleanPars[1]=0.999
760 // 2) Cleaning the tracks
761 // ----------------------
762 // If track's transverse parameter is larger than cleanPars[2]
764 // track's longitudinal parameter is larger than cleanPars[3]
765 // an attempt to remove this track from ESD is made.
767 // The track gets removed if it does not come
768 // from a reconstructed decay
772 Float_t dcaMax=cleanPars[0];
773 Float_t cspMin=cleanPars[1];
775 Int_t nV0s=GetNumberOfV0s();
776 for (Int_t i=nV0s-1; i>=0; i--) {
777 AliESDv0 *v0=GetV0(i);
779 Float_t dca=v0->GetDcaV0Daughters();
780 Float_t csp=v0->GetV0CosineOfPointingAngle();
781 Float_t cspcut=cspMin + dca/dcaMax*(1.-cspMin);
782 if (csp > cspcut) continue;
783 if (RemoveV0(i)) rc=kTRUE;
787 Float_t dmax=cleanPars[2], zmax=cleanPars[3];
789 const AliESDVertex *vertex=GetPrimaryVertexSPD();
790 Bool_t vtxOK=vertex->GetStatus();
792 Int_t nTracks=GetNumberOfTracks();
793 for (Int_t i=nTracks-1; i>=0; i--) {
794 AliESDtrack *track=GetTrack(i);
795 Float_t xy,z; track->GetImpactParameters(xy,z);
796 if ((TMath::Abs(xy) > dmax) || (vtxOK && (TMath::Abs(z) > zmax))) {
797 if (RemoveTrack(i)) rc=kTRUE;
804 Int_t AliESDEvent::AddTrack(const AliESDtrack *t)
807 TClonesArray &ftr = *fTracks;
808 AliESDtrack * track = new(ftr[fTracks->GetEntriesFast()])AliESDtrack(*t);
809 track->SetID(fTracks->GetEntriesFast()-1);
810 return track->GetID();
813 void AliESDEvent::AddMuonTrack(const AliESDMuonTrack *t)
815 TClonesArray &fmu = *fMuonTracks;
816 new(fmu[fMuonTracks->GetEntriesFast()]) AliESDMuonTrack(*t);
819 void AliESDEvent::AddPmdTrack(const AliESDPmdTrack *t)
821 TClonesArray &fpmd = *fPmdTracks;
822 new(fpmd[fPmdTracks->GetEntriesFast()]) AliESDPmdTrack(*t);
825 void AliESDEvent::AddTrdTrack(const AliESDTrdTrack *t)
827 TClonesArray &ftrd = *fTrdTracks;
828 new(ftrd[fTrdTracks->GetEntriesFast()]) AliESDTrdTrack(*t);
834 Int_t AliESDEvent::AddKink(const AliESDkink *c)
837 TClonesArray &fk = *fKinks;
838 AliESDkink * kink = new(fk[fKinks->GetEntriesFast()]) AliESDkink(*c);
839 kink->SetID(fKinks->GetEntriesFast()); // CKB different from the other imps..
840 return fKinks->GetEntriesFast()-1;
844 void AliESDEvent::AddCascade(const AliESDcascade *c)
846 TClonesArray &fc = *fCascades;
847 new(fc[fCascades->GetEntriesFast()]) AliESDcascade(*c);
851 Int_t AliESDEvent::AddCaloCluster(const AliESDCaloCluster *c)
854 TClonesArray &fc = *fCaloClusters;
855 AliESDCaloCluster *clus = new(fc[fCaloClusters->GetEntriesFast()]) AliESDCaloCluster(*c);
856 clus->SetID(fCaloClusters->GetEntriesFast()-1);
857 return fCaloClusters->GetEntriesFast()-1;
861 void AliESDEvent::AddRawDataErrorLog(const AliRawDataErrorLog *log) const {
862 TClonesArray &errlogs = *fErrorLogs;
863 new(errlogs[errlogs.GetEntriesFast()]) AliRawDataErrorLog(*log);
866 void AliESDEvent::SetPrimaryVertexTPC(const AliESDVertex *vertex)
868 // Set the TPC vertex
869 // use already allocated space
871 *fTPCVertex = *vertex;
872 fTPCVertex->SetName(fgkESDListName[kTPCVertex]);
876 void AliESDEvent::SetPrimaryVertexSPD(const AliESDVertex *vertex)
878 // Set the SPD vertex
879 // use already allocated space
881 *fSPDVertex = *vertex;
882 fSPDVertex->SetName(fgkESDListName[kSPDVertex]);
886 void AliESDEvent::SetPrimaryVertex(const AliESDVertex *vertex)
888 // Set the primary vertex
889 // use already allocated space
891 *fPrimaryVertex = *vertex;
892 fPrimaryVertex->SetName(fgkESDListName[kPrimaryVertex]);
896 void AliESDEvent::SetMultiplicity(const AliMultiplicity *mul)
898 // Set the SPD Multiplicity
905 void AliESDEvent::SetFMDData(AliESDFMD * obj)
907 // use already allocated space
913 void AliESDEvent::SetVZEROData(AliESDVZERO * obj)
915 // use already allocated space
920 void AliESDEvent::SetACORDEData(AliESDACORDE * obj)
927 void AliESDEvent::GetESDfriend(AliESDfriend *ev) const
930 // Extracts the complementary info from the ESD
934 Int_t ntrk=GetNumberOfTracks();
936 for (Int_t i=0; i<ntrk; i++) {
937 AliESDtrack *t=GetTrack(i);
938 const AliESDfriendTrack *f=t->GetFriendTrack();
941 t->ReleaseESDfriendTrack();// Not to have two copies of "friendTrack"
945 AliESDfriend *fr = (AliESDfriend*)(const_cast<AliESDEvent*>(this)->FindListObject("AliESDfriend"));
946 if (fr) ev->SetVZEROfriend(fr->GetVZEROfriend());
949 void AliESDEvent::AddObject(TObject* obj)
951 // Add an object to the list of object.
952 // Please be aware that in order to increase performance you should
953 // refrain from using TObjArrays (if possible). Use TClonesArrays, instead.
954 fESDObjects->SetOwner(kTRUE);
955 fESDObjects->AddLast(obj);
959 void AliESDEvent::GetStdContent()
961 // set pointers for standard content
962 // get by name much safer and not a big overhead since not called very often
964 fESDRun = (AliESDRun*)fESDObjects->FindObject(fgkESDListName[kESDRun]);
965 fHeader = (AliESDHeader*)fESDObjects->FindObject(fgkESDListName[kHeader]);
966 fESDZDC = (AliESDZDC*)fESDObjects->FindObject(fgkESDListName[kESDZDC]);
967 fESDFMD = (AliESDFMD*)fESDObjects->FindObject(fgkESDListName[kESDFMD]);
968 fESDVZERO = (AliESDVZERO*)fESDObjects->FindObject(fgkESDListName[kESDVZERO]);
969 fESDTZERO = (AliESDTZERO*)fESDObjects->FindObject(fgkESDListName[kESDTZERO]);
970 fTPCVertex = (AliESDVertex*)fESDObjects->FindObject(fgkESDListName[kTPCVertex]);
971 fSPDVertex = (AliESDVertex*)fESDObjects->FindObject(fgkESDListName[kSPDVertex]);
972 fPrimaryVertex = (AliESDVertex*)fESDObjects->FindObject(fgkESDListName[kPrimaryVertex]);
973 fSPDMult = (AliMultiplicity*)fESDObjects->FindObject(fgkESDListName[kSPDMult]);
974 fPHOSTrigger = (AliESDCaloTrigger*)fESDObjects->FindObject(fgkESDListName[kPHOSTrigger]);
975 fEMCALTrigger = (AliESDCaloTrigger*)fESDObjects->FindObject(fgkESDListName[kEMCALTrigger]);
976 fTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kTracks]);
977 fMuonTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kMuonTracks]);
978 fPmdTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kPmdTracks]);
979 fTrdTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kTrdTracks]);
980 fV0s = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kV0s]);
981 fCascades = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kCascades]);
982 fKinks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kKinks]);
983 fCaloClusters = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kCaloClusters]);
984 fEMCALCells = (AliESDCaloCells*)fESDObjects->FindObject(fgkESDListName[kEMCALCells]);
985 fPHOSCells = (AliESDCaloCells*)fESDObjects->FindObject(fgkESDListName[kPHOSCells]);
986 fErrorLogs = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kErrorLogs]);
987 fESDACORDE = (AliESDACORDE*)fESDObjects->FindObject(fgkESDListName[kESDACORDE]);
991 void AliESDEvent::SetStdNames(){
992 // Set the names of the standard contents
994 if(fESDObjects->GetEntries()>=kESDListN){
995 for(int i = 0;i < fESDObjects->GetEntries() && i<kESDListN;i++){
996 TObject *fObj = fESDObjects->At(i);
997 if(fObj->InheritsFrom("TNamed")){
998 ((TNamed*)fObj)->SetName(fgkESDListName[i]);
1000 else if(fObj->InheritsFrom("TClonesArray")){
1001 ((TClonesArray*)fObj)->SetName(fgkESDListName[i]);
1006 AliWarning("Std Entries missing");
1011 void AliESDEvent::CreateStdContent(Bool_t bUseThisList){
1012 fUseOwnList = bUseThisList;
1016 void AliESDEvent::CreateStdContent()
1018 // create the standard AOD content and set pointers
1020 // create standard objects and add them to the TList of objects
1021 AddObject(new AliESDRun());
1022 AddObject(new AliESDHeader());
1023 AddObject(new AliESDZDC());
1024 AddObject(new AliESDFMD());
1025 AddObject(new AliESDVZERO());
1026 AddObject(new AliESDTZERO());
1027 AddObject(new AliESDVertex());
1028 AddObject(new AliESDVertex());
1029 AddObject(new AliESDVertex());
1030 AddObject(new AliMultiplicity());
1031 AddObject(new AliESDCaloTrigger());
1032 AddObject(new AliESDCaloTrigger());
1033 AddObject(new TClonesArray("AliESDtrack",0));
1034 AddObject(new TClonesArray("AliESDMuonTrack",0));
1035 AddObject(new TClonesArray("AliESDPmdTrack",0));
1036 AddObject(new TClonesArray("AliESDTrdTrack",0));
1037 AddObject(new TClonesArray("AliESDv0",0));
1038 AddObject(new TClonesArray("AliESDcascade",0));
1039 AddObject(new TClonesArray("AliESDkink",0));
1040 AddObject(new TClonesArray("AliESDCaloCluster",0));
1041 AddObject(new AliESDCaloCells());
1042 AddObject(new AliESDCaloCells());
1043 AddObject(new TClonesArray("AliRawDataErrorLog",0));
1044 AddObject(new AliESDACORDE());
1046 // check the order of the indices against enum...
1050 // read back pointers
1054 TObject* AliESDEvent::FindListObject(const char *name){
1056 // Find object with name "name" in the list of branches
1059 return fESDObjects->FindObject(name);
1064 Int_t AliESDEvent::GetPHOSClusters(TRefArray *clusters) const
1066 // fills the provided TRefArray with all found phos clusters
1070 AliESDCaloCluster *cl = 0;
1071 for (Int_t i = 0; i < GetNumberOfCaloClusters(); i++) {
1073 if ( (cl = GetCaloCluster(i)) ) {
1076 AliDebug(1,Form("IsPHOS cluster %d Size: %d \n",i,clusters->GetEntriesFast()));
1080 return clusters->GetEntriesFast();
1083 Int_t AliESDEvent::GetEMCALClusters(TRefArray *clusters) const
1085 // fills the provided TRefArray with all found emcal clusters
1089 AliESDCaloCluster *cl = 0;
1090 for (Int_t i = 0; i < GetNumberOfCaloClusters(); i++) {
1092 if ( (cl = GetCaloCluster(i)) ) {
1095 AliDebug(1,Form("IsEMCAL cluster %d Size: %d \n",i,clusters->GetEntriesFast()));
1099 return clusters->GetEntriesFast();
1102 void AliESDEvent::WriteToTree(TTree* tree) const {
1103 // Book the branches as in TTree::Branch(TCollection*)
1104 // but add a "." at the end of top level branches which are
1105 // not a TClonesArray
1109 TIter next(fESDObjects);
1110 const Int_t kSplitlevel = 99; // default value in TTree::Branch()
1111 const Int_t kBufsize = 32000; // default value in TTree::Branch()
1114 while ((obj = next())) {
1115 branchname.Form("%s", obj->GetName());
1116 if ((kSplitlevel > 1) && !obj->InheritsFrom(TClonesArray::Class())) {
1117 if(!branchname.EndsWith("."))branchname += ".";
1119 if (!tree->FindBranch(branchname)) {
1120 tree->Bronch(branchname, obj->ClassName(), fESDObjects->GetObjectRef(obj),
1121 kBufsize, kSplitlevel - 1);
1127 void AliESDEvent::ReadFromTree(TTree *tree, Option_t* /*opt*/){
1129 // Connect the ESDEvent to a tree
1132 AliWarning("AliESDEvent::ReadFromTree() Zero Pointer to Tree \n");
1136 if(!tree->GetTree())tree->LoadTree(0);
1138 // if we find the "ESD" branch on the tree we do have the old structure
1139 if(tree->GetBranch("ESD")) {
1140 char ** address = (char **)(tree->GetBranch("ESD")->GetAddress());
1141 // do we have the friend branch
1142 TBranch * esdFB = tree->GetBranch("ESDfriend.");
1143 char ** addressF = 0;
1144 if(esdFB)addressF = (char **)(esdFB->GetAddress());
1146 AliInfo("AliESDEvent::ReadFromTree() Reading old Tree");
1147 tree->SetBranchAddress("ESD", &fESDOld);
1149 tree->SetBranchAddress("ESDfriend.",&fESDFriendOld);
1152 AliInfo("AliESDEvent::ReadFromTree() Reading old Tree");
1153 AliInfo("Branch already connected. Using existing branch address.");
1154 fESDOld = (AliESD*) (*address);
1155 // addressF can still be 0, since branch needs to switched on
1156 if(addressF)fESDFriendOld = (AliESDfriend*) (*addressF);
1159 // have already connected the old ESD structure... ?
1160 // reuse also the pointer of the AlliESDEvent
1161 // otherwise create new ones
1162 TList* connectedList = (TList*) (tree->GetUserInfo()->FindObject("ESDObjectsConnectedToTree"));
1165 // If connected use the connected list of objects
1166 if(fESDObjects!= connectedList){
1167 // protect when called twice
1168 fESDObjects->Delete();
1169 fESDObjects = connectedList;
1174 // The pointer to the friend changes when called twice via InitIO
1175 // since AliESDEvent is deleted
1176 TObject* oldf = FindListObject("AliESDfriend");
1179 newf = (TObject*)*addressF;
1181 if(newf!=0&&oldf!=newf){
1182 // remove the old reference
1183 // Should we also delete it? Or is this handled in TTree I/O
1184 // since it is created by the first SetBranchAddress
1185 fESDObjects->Remove(oldf);
1187 fESDObjects->Add(newf);
1194 CreateStdContent(); // create for copy
1195 // if we have the esdfriend add it, so we always can access it via the userinfo
1196 if(fESDFriendOld)AddObject(fESDFriendOld);
1197 // we are not owner of the list objects
1198 // must not delete it
1199 fESDObjects->SetOwner(kFALSE);
1200 fESDObjects->SetName("ESDObjectsConnectedToTree");
1201 tree->GetUserInfo()->Add(fESDObjects);
1209 // Try to find AliESDEvent
1210 AliESDEvent *esdEvent = 0;
1211 esdEvent = (AliESDEvent*)tree->GetTree()->GetUserInfo()->FindObject("AliESDEvent");
1213 // Check if already connected to tree
1214 TList* connectedList = (TList*) (tree->GetUserInfo()->FindObject("ESDObjectsConnectedToTree"));
1215 if (connectedList) {
1216 // If connected use the connected list if objects
1217 fESDObjects->Delete();
1218 fESDObjects = connectedList;
1225 // prevent a memory leak when reading back the TList
1230 // create a new TList from the UserInfo TList...
1231 // copy constructor does not work...
1232 fESDObjects = (TList*)(esdEvent->GetList()->Clone());
1233 fESDObjects->SetOwner(kFALSE);
1235 else if ( fESDObjects->GetEntries()==0){
1236 // at least create the std content if we want to read to our list
1241 // we only need new things in the list if we do no already have it..
1242 // TODO just add new entries
1244 if(fESDObjects->GetEntries()<kESDListN){
1245 AliWarning(Form("AliESDEvent::ReadFromTree() TList contains less than the standard contents %d < %d \n",
1246 fESDObjects->GetEntries(),kESDListN));
1248 // set the branch addresses
1249 TIter next(fESDObjects);
1251 while((el=(TNamed*)next())){
1252 TString bname(el->GetName());
1253 if(bname.CompareTo("AliESDfriend")==0)
1255 // AliESDfriend does not have a name ...
1256 tree->SetBranchAddress("ESDfriend.",fESDObjects->GetObjectRef(el));
1259 // check if branch exists under this Name
1260 TBranch *br = tree->GetBranch(bname.Data());
1262 tree->SetBranchAddress(bname.Data(),fESDObjects->GetObjectRef(el));
1265 br = tree->GetBranch(Form("%s.",bname.Data()));
1267 tree->SetBranchAddress(Form("%s.",bname.Data()),fESDObjects->GetObjectRef(el));
1270 AliWarning(Form("AliESDEvent::ReadFromTree() No Branch found with Name %s or %s.",bname.Data(),bname.Data()));
1277 // when reading back we are not owner of the list
1278 // must not delete it
1279 fESDObjects->SetOwner(kFALSE);
1280 fESDObjects->SetName("ESDObjectsConnectedToTree");
1281 // we are not owner of the list objects
1282 // must not delete it
1283 tree->GetUserInfo()->Add(fESDObjects);
1287 // we can't get the list from the user data, create standard content
1288 // and set it by hand (no ESDfriend at the moment
1290 TIter next(fESDObjects);
1292 while((el=(TNamed*)next())){
1293 TString bname(el->GetName());
1294 TBranch *br = tree->GetBranch(bname.Data());
1296 tree->SetBranchAddress(bname.Data(),fESDObjects->GetObjectRef(el));
1299 br = tree->GetBranch(Form("%s.",bname.Data()));
1301 tree->SetBranchAddress(Form("%s.",bname.Data()),fESDObjects->GetObjectRef(el));
1306 // when reading back we are not owner of the list
1307 // must not delete it
1308 fESDObjects->SetOwner(kFALSE);
1313 void AliESDEvent::CopyFromOldESD()
1315 // Method which copies over everthing from the old esd structure to the
1320 SetRunNumber(fESDOld->GetRunNumber());
1321 SetPeriodNumber(fESDOld->GetPeriodNumber());
1322 SetMagneticField(fESDOld->GetMagneticField());
1324 // leave out diamond ...
1325 // SetDiamond(const AliESDVertex *vertex) { fESDRun->SetDiamond(vertex);}
1328 SetTriggerMask(fESDOld->GetTriggerMask());
1329 SetOrbitNumber(fESDOld->GetOrbitNumber());
1330 SetTimeStamp(fESDOld->GetTimeStamp());
1331 SetEventType(fESDOld->GetEventType());
1332 SetEventNumberInFile(fESDOld->GetEventNumberInFile());
1333 SetBunchCrossNumber(fESDOld->GetBunchCrossNumber());
1334 SetTriggerCluster(fESDOld->GetTriggerCluster());
1338 SetZDC(fESDOld->GetZDCN1Energy(),
1339 fESDOld->GetZDCP1Energy(),
1340 fESDOld->GetZDCEMEnergy(),
1342 fESDOld->GetZDCN2Energy(),
1343 fESDOld->GetZDCP2Energy(),
1344 fESDOld->GetZDCParticipants(),
1349 if(fESDOld->GetFMDData())SetFMDData(fESDOld->GetFMDData());
1353 SetT0zVertex(fESDOld->GetT0zVertex());
1354 SetT0(fESDOld->GetT0());
1358 if (fESDOld->GetVZEROData()) SetVZEROData(fESDOld->GetVZEROData());
1360 if(fESDOld->GetVertex())SetPrimaryVertexSPD(fESDOld->GetVertex());
1362 if(fESDOld->GetPrimaryVertex())SetPrimaryVertex(fESDOld->GetPrimaryVertex());
1364 if(fESDOld->GetMultiplicity())SetMultiplicity(fESDOld->GetMultiplicity());
1366 for(int i = 0;i<fESDOld->GetNumberOfTracks();i++){
1367 AddTrack(fESDOld->GetTrack(i));
1370 for(int i = 0;i<fESDOld->GetNumberOfMuonTracks();i++){
1371 AddMuonTrack(fESDOld->GetMuonTrack(i));
1374 for(int i = 0;i<fESDOld->GetNumberOfPmdTracks();i++){
1375 AddPmdTrack(fESDOld->GetPmdTrack(i));
1378 for(int i = 0;i<fESDOld->GetNumberOfTrdTracks();i++){
1379 AddTrdTrack(fESDOld->GetTrdTrack(i));
1382 for(int i = 0;i<fESDOld->GetNumberOfV0s();i++){
1383 AddV0(fESDOld->GetV0(i));
1386 for(int i = 0;i<fESDOld->GetNumberOfCascades();i++){
1387 AddCascade(fESDOld->GetCascade(i));
1390 for(int i = 0;i<fESDOld->GetNumberOfKinks();i++){
1391 AddKink(fESDOld->GetKink(i));
1395 for(int i = 0;i<fESDOld->GetNumberOfCaloClusters();i++){
1396 AddCaloCluster(fESDOld->GetCaloCluster(i));