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),
127 fFirstEMCALCluster(-1),
129 fFirstPHOSCluster(-1)
132 //______________________________________________________________________________
133 AliESDEvent::AliESDEvent(const AliESDEvent& esd):
135 fESDObjects(new TList()),
136 fESDRun(new AliESDRun(*esd.fESDRun)),
137 fHeader(new AliESDHeader(*esd.fHeader)),
138 fESDZDC(new AliESDZDC(*esd.fESDZDC)),
139 fESDFMD(new AliESDFMD(*esd.fESDFMD)),
140 fESDVZERO(new AliESDVZERO(*esd.fESDVZERO)),
141 fESDTZERO(new AliESDTZERO(*esd.fESDTZERO)),
142 fTPCVertex(new AliESDVertex(*esd.fTPCVertex)),
143 fSPDVertex(new AliESDVertex(*esd.fSPDVertex)),
144 fPrimaryVertex(new AliESDVertex(*esd.fPrimaryVertex)),
145 fSPDMult(new AliMultiplicity(*esd.fSPDMult)),
146 fPHOSTrigger(new AliESDCaloTrigger(*esd.fPHOSTrigger)),
147 fEMCALTrigger(new AliESDCaloTrigger(*esd.fEMCALTrigger)),
148 fESDACORDE(new AliESDACORDE(*esd.fESDACORDE)),
149 fTracks(new TClonesArray(*esd.fTracks)),
150 fMuonTracks(new TClonesArray(*esd.fMuonTracks)),
151 fPmdTracks(new TClonesArray(*esd.fPmdTracks)),
152 fTrdTracks(new TClonesArray(*esd.fTrdTracks)),
153 fV0s(new TClonesArray(*esd.fV0s)),
154 fCascades(new TClonesArray(*esd.fCascades)),
155 fKinks(new TClonesArray(*esd.fKinks)),
156 fCaloClusters(new TClonesArray(*esd.fCaloClusters)),
157 fEMCALCells(new AliESDCaloCells(*esd.fEMCALCells)),
158 fPHOSCells(new AliESDCaloCells(*esd.fPHOSCells)),
159 fErrorLogs(new TClonesArray(*esd.fErrorLogs)),
160 fESDOld(new AliESD(*esd.fESDOld)),
161 fESDFriendOld(new AliESDfriend(*esd.fESDFriendOld)),
162 fConnected(esd.fConnected),
163 fEMCALClusters(esd.fEMCALClusters),
164 fFirstEMCALCluster(esd.fFirstEMCALCluster),
165 fPHOSClusters(esd.fPHOSClusters),
166 fFirstPHOSCluster(esd.fFirstPHOSCluster)
169 // CKB init in the constructor list and only add here ...
174 AddObject(fESDVZERO);
175 AddObject(fESDTZERO);
176 AddObject(fTPCVertex);
177 AddObject(fSPDVertex);
178 AddObject(fPrimaryVertex);
180 AddObject(fPHOSTrigger);
181 AddObject(fEMCALTrigger);
183 AddObject(fMuonTracks);
184 AddObject(fPmdTracks);
185 AddObject(fTrdTracks);
187 AddObject(fCascades);
189 AddObject(fCaloClusters);
190 AddObject(fEMCALCells);
191 AddObject(fPHOSCells);
192 AddObject(fErrorLogs);
193 AddObject(fESDACORDE);
199 //______________________________________________________________________________
200 AliESDEvent & AliESDEvent::operator=(const AliESDEvent& source) {
202 // Assignment operator
204 if(&source == this) return *this;
205 AliVEvent::operator=(source);
207 fESDRun = new AliESDRun(*source.fESDRun);
208 fHeader = new AliESDHeader(*source.fHeader);
209 fESDZDC = new AliESDZDC(*source.fESDZDC);
210 fESDFMD = new AliESDFMD(*source.fESDFMD);
211 fESDVZERO = new AliESDVZERO(*source.fESDVZERO);
212 fESDTZERO = new AliESDTZERO(*source.fESDTZERO);
213 fTPCVertex = new AliESDVertex(*source.fTPCVertex);
214 fSPDVertex = new AliESDVertex(*source.fSPDVertex);
215 fPrimaryVertex = new AliESDVertex(*source.fPrimaryVertex);
216 fSPDMult = new AliMultiplicity(*source.fSPDMult);
217 fPHOSTrigger = new AliESDCaloTrigger(*source.fPHOSTrigger);
218 fEMCALTrigger = new AliESDCaloTrigger(*source.fEMCALTrigger);
219 fTracks = new TClonesArray(*source.fTracks);
220 fMuonTracks = new TClonesArray(*source.fMuonTracks);
221 fPmdTracks = new TClonesArray(*source.fPmdTracks);
222 fTrdTracks = new TClonesArray(*source.fTrdTracks);
223 fV0s = new TClonesArray(*source.fV0s);
224 fCascades = new TClonesArray(*source.fCascades);
225 fKinks = new TClonesArray(*source.fKinks);
226 fCaloClusters = new TClonesArray(*source.fCaloClusters);
227 fEMCALCells = new AliESDCaloCells(*source.fEMCALCells);
228 fPHOSCells = new AliESDCaloCells(*source.fPHOSCells);
229 fErrorLogs = new TClonesArray(*source.fErrorLogs);
230 fESDACORDE = new AliESDACORDE(*source.fESDACORDE);
231 fESDOld = new AliESD(*source.fESDOld);
232 fESDFriendOld = new AliESDfriend(*source.fESDFriendOld);
234 // or AddObject( fESDZDC = new AliESDZDC(*source.fESDZDC));
236 fESDObjects = new TList();
241 AddObject(fESDVZERO);
242 AddObject(fESDTZERO);
243 AddObject(fTPCVertex);
244 AddObject(fSPDVertex);
245 AddObject(fPrimaryVertex);
247 AddObject(fPHOSTrigger);
248 AddObject(fEMCALTrigger);
250 AddObject(fMuonTracks);
251 AddObject(fPmdTracks);
252 AddObject(fTrdTracks);
254 AddObject(fCascades);
256 AddObject(fCaloClusters);
257 AddObject(fEMCALCells);
258 AddObject(fPHOSCells);
259 AddObject(fErrorLogs);
260 AddObject(fESDACORDE);
262 fConnected = source.fConnected;
263 fEMCALClusters = source.fEMCALClusters;
264 fFirstEMCALCluster = source.fFirstEMCALCluster;
265 fPHOSClusters = source.fPHOSClusters;
266 fFirstPHOSCluster = source.fFirstPHOSCluster;
275 //______________________________________________________________________________
276 AliESDEvent::~AliESDEvent()
279 // Standard destructor
282 // everthing on the list gets deleted automatically
285 if(fESDObjects&&!fConnected)
294 //______________________________________________________________________________
295 void AliESDEvent::Reset()
299 // Reset the standard contents
301 if(fESDOld)fESDOld->Reset();
302 // reset for the friends...
304 fESDFriendOld->~AliESDfriend();
305 new (fESDFriendOld) AliESDfriend();
307 // for new data we have to fetch the Pointer from the list
308 AliESDfriend *fr = (AliESDfriend*)FindListObject("AliESDfriend");
310 // delete the content
312 // make a new valid ESDfriend at the same place
313 new (fr) AliESDfriend();
316 // call reset for user supplied data?
319 void AliESDEvent::ResetStdContent()
321 // Reset the standard contents
322 if(fESDRun) fESDRun->Reset();
323 if(fHeader) fHeader->Reset();
324 if(fESDZDC) fESDZDC->Reset();
326 fESDFMD->~AliESDFMD();
327 new (fESDFMD) AliESDFMD();
330 // reset by callin d'to /c'tor keep the pointer
331 fESDVZERO->~AliESDVZERO();
332 new (fESDVZERO) AliESDVZERO();
335 fESDACORDE->~AliESDACORDE();
336 new (fESDACORDE) AliESDACORDE();
338 if(fESDTZERO) fESDTZERO->Reset();
339 // CKB no clear/reset implemented
341 fTPCVertex->~AliESDVertex();
342 new (fTPCVertex) AliESDVertex();
343 fTPCVertex->SetName(fgkESDListName[kTPCVertex]);
346 fSPDVertex->~AliESDVertex();
347 new (fSPDVertex) AliESDVertex();
348 fSPDVertex->SetName(fgkESDListName[kSPDVertex]);
351 fPrimaryVertex->~AliESDVertex();
352 new (fPrimaryVertex) AliESDVertex();
353 fPrimaryVertex->SetName(fgkESDListName[kPrimaryVertex]);
356 fSPDMult->~AliMultiplicity();
357 new (fSPDMult) AliMultiplicity();
359 if(fPHOSTrigger)fPHOSTrigger->Reset();
360 if(fEMCALTrigger)fEMCALTrigger->Reset();
361 if(fTracks)fTracks->Delete();
362 if(fMuonTracks)fMuonTracks->Delete();
363 if(fPmdTracks)fPmdTracks->Delete();
364 if(fTrdTracks)fTrdTracks->Delete();
365 if(fV0s)fV0s->Delete();
366 if(fCascades)fCascades->Delete();
367 if(fKinks)fKinks->Delete();
368 if(fCaloClusters)fCaloClusters->Delete();
369 if(fPHOSCells)fPHOSCells->DeleteContainer();
370 if(fEMCALCells)fEMCALCells->DeleteContainer();
371 if(fErrorLogs) fErrorLogs->Delete();
373 // don't reset fconnected fConnected and the list
376 fFirstEMCALCluster=-1;
378 fFirstPHOSCluster=-1;
382 Int_t AliESDEvent::AddV0(const AliESDv0 *v) {
386 TClonesArray &fv = *fV0s;
387 Int_t idx=fV0s->GetEntriesFast();
388 new(fv[idx]) AliESDv0(*v);
392 //______________________________________________________________________________
393 void AliESDEvent::Print(Option_t *) const
396 // Print header information of the event
398 printf("ESD run information\n");
399 printf("Event # in file %d Bunch crossing # %d Orbit # %d Period # %d Run # %d Trigger %lld Magnetic field %f \n",
400 GetEventNumberInFile(),
401 GetBunchCrossNumber(),
406 GetMagneticField() );
407 printf("Vertex: (%.4f +- %.4f, %.4f +- %.4f, %.4f +- %.4f) cm\n",
408 fPrimaryVertex->GetXv(), fPrimaryVertex->GetXRes(),
409 fPrimaryVertex->GetYv(), fPrimaryVertex->GetYRes(),
410 fPrimaryVertex->GetZv(), fPrimaryVertex->GetZRes());
411 printf("Mean vertex in RUN: X=%.4f Y=%.4f cm\n",
412 GetDiamondX(),GetDiamondY());
413 printf("SPD Multiplicity. Number of tracklets %d \n",
414 fSPDMult->GetNumberOfTracklets());
415 printf("Number of tracks: \n");
416 printf(" charged %d\n", GetNumberOfTracks());
417 printf(" muon %d\n", GetNumberOfMuonTracks());
418 printf(" pmd %d\n", GetNumberOfPmdTracks());
419 printf(" trd %d\n", GetNumberOfTrdTracks());
420 printf(" v0 %d\n", GetNumberOfV0s());
421 printf(" cascades %d\n", GetNumberOfCascades());
422 printf(" kinks %d\n", GetNumberOfKinks());
423 if(fPHOSCells)printf(" PHOSCells %d\n", fPHOSCells->GetNumberOfCells());
424 else printf(" PHOSCells not in the Event\n");
425 if(fEMCALCells)printf(" EMCALCells %d\n", fEMCALCells->GetNumberOfCells());
426 else printf(" EMCALCells not in the Event\n");
427 printf(" CaloClusters %d\n", GetNumberOfCaloClusters());
428 printf(" phos %d\n", GetNumberOfPHOSClusters());
429 printf(" emcal %d\n", GetNumberOfEMCALClusters());
430 printf(" FMD %s\n", (fESDFMD ? "yes" : "no"));
431 printf(" VZERO %s\n", (fESDVZERO ? "yes" : "no"));
436 void AliESDEvent::SetESDfriend(const AliESDfriend *ev) const {
438 // Attaches the complementary info to the ESD
442 // to be sure that we set the tracks also
443 // in case of old esds
444 // if(fESDOld)CopyFromOldESD();
446 Int_t ntrk=ev->GetNumberOfTracks();
448 for (Int_t i=0; i<ntrk; i++) {
449 const AliESDfriendTrack *f=ev->GetTrack(i);
450 GetTrack(i)->SetFriendTrack(f);
454 Bool_t AliESDEvent::RemoveKink(Int_t rm) const {
455 // ---------------------------------------------------------
456 // Remove a kink candidate and references to it from ESD,
457 // if this candidate does not come from a reconstructed decay
458 // Not yet implemented...
459 // ---------------------------------------------------------
460 Int_t last=GetNumberOfKinks()-1;
461 if ((rm<0)||(rm>last)) return kFALSE;
466 Bool_t AliESDEvent::RemoveV0(Int_t rm) const {
467 // ---------------------------------------------------------
468 // Remove a V0 candidate and references to it from ESD,
469 // if this candidate does not come from a reconstructed decay
470 // ---------------------------------------------------------
471 Int_t last=GetNumberOfV0s()-1;
472 if ((rm<0)||(rm>last)) return kFALSE;
474 AliESDv0 *v0=GetV0(rm);
475 Int_t idxP=v0->GetPindex(), idxN=v0->GetNindex();
478 Int_t lastIdxP=v0->GetPindex(), lastIdxN=v0->GetNindex();
482 // Check if this V0 comes from a reconstructed decay
483 Int_t ncs=GetNumberOfCascades();
484 for (Int_t n=0; n<ncs; n++) {
485 AliESDcascade *cs=GetCascade(n);
487 Int_t csIdxP=cs->GetPindex();
488 Int_t csIdxN=cs->GetNindex();
491 if (idxN==csIdxN) return kFALSE;
493 if (csIdxP==lastIdxP)
494 if (csIdxN==lastIdxN) used++;
497 //Replace the removed V0 with the last V0
498 TClonesArray &a=*fV0s;
499 delete a.RemoveAt(rm);
501 if (rm==last) return kTRUE;
503 //v0 is pointing to the last V0 candidate...
504 new (a[rm]) AliESDv0(*v0);
505 delete a.RemoveAt(last);
507 if (!used) return kTRUE;
510 // Remap the indices of the daughters of reconstructed decays
511 for (Int_t n=0; n<ncs; n++) {
512 AliESDcascade *cs=GetCascade(n);
515 Int_t csIdxP=cs->GetPindex();
516 Int_t csIdxN=cs->GetNindex();
518 if (csIdxP==lastIdxP)
519 if (csIdxN==lastIdxN) {
520 cs->AliESDv0::SetIndex(1,idxP);
521 cs->AliESDv0::SetIndex(0,idxN);
523 if (!used) return kTRUE;
530 Bool_t AliESDEvent::RemoveTrack(Int_t rm) const {
531 // ---------------------------------------------------------
532 // Remove a track and references to it from ESD,
533 // if this track does not come from a reconstructed decay
534 // ---------------------------------------------------------
535 Int_t last=GetNumberOfTracks()-1;
536 if ((rm<0)||(rm>last)) return kFALSE;
540 // Check if this track comes from the reconstructed primary vertices
541 if (fTPCVertex && fTPCVertex->GetStatus()) {
542 UShort_t *primIdx=fTPCVertex->GetIndices();
543 Int_t n=fTPCVertex->GetNIndices();
545 Int_t idx=Int_t(primIdx[n]);
546 if (rm==idx) return kFALSE;
547 if (idx==last) used++;
550 if (fPrimaryVertex && fPrimaryVertex->GetStatus()) {
551 UShort_t *primIdx=fPrimaryVertex->GetIndices();
552 Int_t n=fPrimaryVertex->GetNIndices();
554 Int_t idx=Int_t(primIdx[n]);
555 if (rm==idx) return kFALSE;
556 if (idx==last) used++;
560 // Check if this track comes from a reconstructed decay
561 Int_t nv0=GetNumberOfV0s();
562 for (Int_t n=0; n<nv0; n++) {
563 AliESDv0 *v0=GetV0(n);
565 Int_t idx=v0->GetNindex();
566 if (rm==idx) return kFALSE;
567 if (idx==last) used++;
570 if (rm==idx) return kFALSE;
571 if (idx==last) used++;
574 Int_t ncs=GetNumberOfCascades();
575 for (Int_t n=0; n<ncs; n++) {
576 AliESDcascade *cs=GetCascade(n);
578 Int_t idx=cs->GetIndex();
579 if (rm==idx) return kFALSE;
580 if (idx==last) used++;
583 Int_t nkn=GetNumberOfKinks();
584 for (Int_t n=0; n<nkn; n++) {
585 AliESDkink *kn=GetKink(n);
587 Int_t idx=kn->GetIndex(0);
588 if (rm==idx) return kFALSE;
589 if (idx==last) used++;
592 if (rm==idx) return kFALSE;
593 if (idx==last) used++;
597 //Replace the removed track with the last track
598 TClonesArray &a=*fTracks;
599 delete a.RemoveAt(rm);
601 if (rm==last) return kTRUE;
603 AliESDtrack *t=GetTrack(last);
605 new (a[rm]) AliESDtrack(*t);
606 delete a.RemoveAt(last);
609 if (!used) return kTRUE;
612 // Remap the indices of the tracks used for the primary vertex reconstruction
613 if (fTPCVertex && fTPCVertex->GetStatus()) {
614 UShort_t *primIdx=fTPCVertex->GetIndices();
615 Int_t n=fTPCVertex->GetNIndices();
617 Int_t idx=Int_t(primIdx[n]);
619 primIdx[n]=Short_t(rm);
621 if (!used) return kTRUE;
625 if (fPrimaryVertex && fPrimaryVertex->GetStatus()) {
626 UShort_t *primIdx=fPrimaryVertex->GetIndices();
627 Int_t n=fPrimaryVertex->GetNIndices();
629 Int_t idx=Int_t(primIdx[n]);
631 primIdx[n]=Short_t(rm);
633 if (!used) return kTRUE;
638 // Remap the indices of the daughters of reconstructed decays
639 for (Int_t n=0; n<nv0; n++) {
640 AliESDv0 *v0=GetV0(n);
641 if (v0->GetIndex(0)==last) {
644 if (!used) return kTRUE;
646 if (v0->GetIndex(1)==last) {
649 if (!used) return kTRUE;
653 for (Int_t n=0; n<ncs; n++) {
654 AliESDcascade *cs=GetCascade(n);
655 if (cs->GetIndex()==last) {
658 if (!used) return kTRUE;
662 for (Int_t n=0; n<nkn; n++) {
663 AliESDkink *kn=GetKink(n);
664 if (kn->GetIndex(0)==last) {
667 if (!used) return kTRUE;
669 if (kn->GetIndex(1)==last) {
672 if (!used) return kTRUE;
680 Bool_t AliESDEvent::Clean(Float_t *cleanPars) {
682 // Remove the data which are not needed for the physics analysis.
684 // 1) Cleaning the V0 candidates
685 // ---------------------------
686 // If the cosine of the V0 pointing angle "csp" and
687 // the DCA between the daughter tracks "dca" does not satisfy
690 // csp > cleanPars[1] + dca/cleanPars[0]*(1.- cleanPars[1])
692 // an attempt to remove this V0 candidate from ESD is made.
694 // The V0 candidate gets removed if it does not belong to any
695 // recosntructed cascade decay
697 // 12.11.2007, optimal values: cleanPars[0]=0.5, cleanPars[1]=0.999
699 // 2) Cleaning the tracks
700 // ----------------------
701 // If track's transverse parameter is larger than cleanPars[2]
703 // track's longitudinal parameter is larger than cleanPars[3]
704 // an attempt to remove this track from ESD is made.
706 // The track gets removed if it does not come
707 // from a reconstructed decay
711 Float_t dcaMax=cleanPars[0];
712 Float_t cspMin=cleanPars[1];
714 Int_t nV0s=GetNumberOfV0s();
715 for (Int_t i=nV0s-1; i>=0; i--) {
716 AliESDv0 *v0=GetV0(i);
718 Float_t dca=v0->GetDcaV0Daughters();
719 Float_t csp=v0->GetV0CosineOfPointingAngle();
720 Float_t cspcut=cspMin + dca/dcaMax*(1.-cspMin);
721 if (csp > cspcut) continue;
722 if (RemoveV0(i)) rc=kTRUE;
726 Float_t dmax=cleanPars[2], zmax=cleanPars[3];
728 const AliESDVertex *vertex=GetPrimaryVertexSPD();
729 Bool_t vtxOK=vertex->GetStatus();
731 Int_t nTracks=GetNumberOfTracks();
732 for (Int_t i=nTracks-1; i>=0; i--) {
733 AliESDtrack *track=GetTrack(i);
734 Float_t xy,z; track->GetImpactParameters(xy,z);
735 if ((TMath::Abs(xy) > dmax) || (vtxOK && (TMath::Abs(z) > zmax))) {
736 if (RemoveTrack(i)) rc=kTRUE;
743 Int_t AliESDEvent::AddTrack(const AliESDtrack *t)
746 TClonesArray &ftr = *fTracks;
747 AliESDtrack * track = new(ftr[fTracks->GetEntriesFast()])AliESDtrack(*t);
748 track->SetID(fTracks->GetEntriesFast()-1);
749 return track->GetID();
752 void AliESDEvent::AddMuonTrack(const AliESDMuonTrack *t)
754 TClonesArray &fmu = *fMuonTracks;
755 new(fmu[fMuonTracks->GetEntriesFast()]) AliESDMuonTrack(*t);
758 void AliESDEvent::AddPmdTrack(const AliESDPmdTrack *t)
760 TClonesArray &fpmd = *fPmdTracks;
761 new(fpmd[fPmdTracks->GetEntriesFast()]) AliESDPmdTrack(*t);
764 void AliESDEvent::AddTrdTrack(const AliESDTrdTrack *t)
766 TClonesArray &ftrd = *fTrdTracks;
767 new(ftrd[fTrdTracks->GetEntriesFast()]) AliESDTrdTrack(*t);
773 Int_t AliESDEvent::AddKink(const AliESDkink *c)
776 TClonesArray &fk = *fKinks;
777 AliESDkink * kink = new(fk[fKinks->GetEntriesFast()]) AliESDkink(*c);
778 kink->SetID(fKinks->GetEntriesFast()); // CKB different from the other imps..
779 return fKinks->GetEntriesFast()-1;
783 void AliESDEvent::AddCascade(const AliESDcascade *c)
785 TClonesArray &fc = *fCascades;
786 new(fc[fCascades->GetEntriesFast()]) AliESDcascade(*c);
790 Int_t AliESDEvent::AddCaloCluster(const AliESDCaloCluster *c)
793 TClonesArray &fc = *fCaloClusters;
794 AliESDCaloCluster *clus = new(fc[fCaloClusters->GetEntriesFast()]) AliESDCaloCluster(*c);
795 clus->SetID(fCaloClusters->GetEntriesFast()-1);
796 return fCaloClusters->GetEntriesFast()-1;
800 void AliESDEvent::AddRawDataErrorLog(const AliRawDataErrorLog *log) const {
801 TClonesArray &errlogs = *fErrorLogs;
802 new(errlogs[errlogs.GetEntriesFast()]) AliRawDataErrorLog(*log);
805 void AliESDEvent::SetPrimaryVertexTPC(const AliESDVertex *vertex)
807 // Set the TPC vertex
808 // use already allocated space
810 *fTPCVertex = *vertex;
811 fTPCVertex->SetName(fgkESDListName[kTPCVertex]);
815 void AliESDEvent::SetPrimaryVertexSPD(const AliESDVertex *vertex)
817 // Set the SPD vertex
818 // use already allocated space
820 *fSPDVertex = *vertex;
821 fSPDVertex->SetName(fgkESDListName[kSPDVertex]);
825 void AliESDEvent::SetPrimaryVertex(const AliESDVertex *vertex)
827 // Set the primary vertex
828 // use already allocated space
830 *fPrimaryVertex = *vertex;
831 fPrimaryVertex->SetName(fgkESDListName[kPrimaryVertex]);
835 void AliESDEvent::SetMultiplicity(const AliMultiplicity *mul)
837 // Set the SPD Multiplicity
844 void AliESDEvent::SetFMDData(AliESDFMD * obj)
846 // use already allocated space
852 void AliESDEvent::SetVZEROData(AliESDVZERO * obj)
854 // use already allocated space
859 void AliESDEvent::SetACORDEData(AliESDACORDE * obj)
866 void AliESDEvent::GetESDfriend(AliESDfriend *ev) const
869 // Extracts the complementary info from the ESD
873 Int_t ntrk=GetNumberOfTracks();
875 for (Int_t i=0; i<ntrk; i++) {
876 AliESDtrack *t=GetTrack(i);
877 const AliESDfriendTrack *f=t->GetFriendTrack();
880 t->ReleaseESDfriendTrack();// Not to have two copies of "friendTrack"
886 void AliESDEvent::AddObject(TObject* obj)
888 // Add an object to the list of object.
889 // Please be aware that in order to increase performance you should
890 // refrain from using TObjArrays (if possible). Use TClonesArrays, instead.
891 fESDObjects->SetOwner(kTRUE);
892 fESDObjects->AddLast(obj);
896 void AliESDEvent::GetStdContent()
898 // set pointers for standard content
899 // get by name much safer and not a big overhead since not called very often
901 fESDRun = (AliESDRun*)fESDObjects->FindObject(fgkESDListName[kESDRun]);
902 fHeader = (AliESDHeader*)fESDObjects->FindObject(fgkESDListName[kHeader]);
903 fESDZDC = (AliESDZDC*)fESDObjects->FindObject(fgkESDListName[kESDZDC]);
904 fESDFMD = (AliESDFMD*)fESDObjects->FindObject(fgkESDListName[kESDFMD]);
905 fESDVZERO = (AliESDVZERO*)fESDObjects->FindObject(fgkESDListName[kESDVZERO]);
906 fESDTZERO = (AliESDTZERO*)fESDObjects->FindObject(fgkESDListName[kESDTZERO]);
907 fTPCVertex = (AliESDVertex*)fESDObjects->FindObject(fgkESDListName[kTPCVertex]);
908 fSPDVertex = (AliESDVertex*)fESDObjects->FindObject(fgkESDListName[kSPDVertex]);
909 fPrimaryVertex = (AliESDVertex*)fESDObjects->FindObject(fgkESDListName[kPrimaryVertex]);
910 fSPDMult = (AliMultiplicity*)fESDObjects->FindObject(fgkESDListName[kSPDMult]);
911 fPHOSTrigger = (AliESDCaloTrigger*)fESDObjects->FindObject(fgkESDListName[kPHOSTrigger]);
912 fEMCALTrigger = (AliESDCaloTrigger*)fESDObjects->FindObject(fgkESDListName[kEMCALTrigger]);
913 fTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kTracks]);
914 fMuonTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kMuonTracks]);
915 fPmdTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kPmdTracks]);
916 fTrdTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kTrdTracks]);
917 fV0s = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kV0s]);
918 fCascades = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kCascades]);
919 fKinks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kKinks]);
920 fCaloClusters = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kCaloClusters]);
921 fEMCALCells = (AliESDCaloCells*)fESDObjects->FindObject(fgkESDListName[kEMCALCells]);
922 fPHOSCells = (AliESDCaloCells*)fESDObjects->FindObject(fgkESDListName[kPHOSCells]);
923 fErrorLogs = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kErrorLogs]);
924 fESDACORDE = (AliESDACORDE*)fESDObjects->FindObject(fgkESDListName[kESDACORDE]);
928 void AliESDEvent::SetStdNames(){
929 // Set the names of the standard contents
931 if(fESDObjects->GetEntries()==kESDListN){
932 for(int i = 0;i < fESDObjects->GetEntries();i++){
933 TObject *fObj = fESDObjects->At(i);
934 if(fObj->InheritsFrom("TNamed")){
935 ((TNamed*)fObj)->SetName(fgkESDListName[i]);
937 else if(fObj->InheritsFrom("TClonesArray")){
938 ((TClonesArray*)fObj)->SetName(fgkESDListName[i]);
943 printf("%s:%d SetStdNames() Wrong number of Std Entries \n",(char*)__FILE__,__LINE__);
947 void AliESDEvent::CreateStdContent()
949 // create the standard AOD content and set pointers
951 // create standard objects and add them to the TList of objects
952 AddObject(new AliESDRun());
953 AddObject(new AliESDHeader());
954 AddObject(new AliESDZDC());
955 AddObject(new AliESDFMD());
956 AddObject(new AliESDVZERO());
957 AddObject(new AliESDTZERO());
958 AddObject(new AliESDVertex());
959 AddObject(new AliESDVertex());
960 AddObject(new AliESDVertex());
961 AddObject(new AliMultiplicity());
962 AddObject(new AliESDCaloTrigger());
963 AddObject(new AliESDCaloTrigger());
964 AddObject(new TClonesArray("AliESDtrack",0));
965 AddObject(new TClonesArray("AliESDMuonTrack",0));
966 AddObject(new TClonesArray("AliESDPmdTrack",0));
967 AddObject(new TClonesArray("AliESDTrdTrack",0));
968 AddObject(new TClonesArray("AliESDv0",0));
969 AddObject(new TClonesArray("AliESDcascade",0));
970 AddObject(new TClonesArray("AliESDkink",0));
971 AddObject(new TClonesArray("AliESDCaloCluster",0));
972 AddObject(new AliESDCaloCells());
973 AddObject(new AliESDCaloCells());
974 AddObject(new TClonesArray("AliRawDataErrorLog",0));
975 AddObject(new AliESDACORDE());
977 // check the order of the indices against enum...
981 // read back pointers
985 TObject* AliESDEvent::FindListObject(const char *name){
987 // Find object with name "name" in the list of branches
990 return fESDObjects->FindObject(name);
995 Int_t AliESDEvent::GetPHOSClusters(TRefArray *clusters) const
997 // fills the provided TRefArray with all found phos clusters
1001 AliESDCaloCluster *cl = 0;
1002 for (Int_t i = 0; i < GetNumberOfCaloClusters(); i++) {
1004 if ( (cl = GetCaloCluster(i)) ) {
1007 AliDebug(1,Form("IsPHOS cluster %d Size: %d \n",i,clusters->GetEntriesFast()));
1011 return clusters->GetEntriesFast();
1014 Int_t AliESDEvent::GetEMCALClusters(TRefArray *clusters) const
1016 // fills the provided TRefArray with all found emcal clusters
1020 AliESDCaloCluster *cl = 0;
1021 for (Int_t i = 0; i < GetNumberOfCaloClusters(); i++) {
1023 if ( (cl = GetCaloCluster(i)) ) {
1026 AliDebug(1,Form("IsEMCAL cluster %d Size: %d \n",i,clusters->GetEntriesFast()));
1030 return clusters->GetEntriesFast();
1033 const void AliESDEvent::WriteToTree(TTree* tree) const {
1034 // Book the branches as in TTree::Branch(TCollection*)
1035 // but add a "." at the end of top level branches which are
1036 // not a TClonesArray
1040 TIter next(fESDObjects);
1041 const Int_t kSplitlevel = 99; // default value in TTree::Branch()
1042 const Int_t kBufsize = 32000; // default value in TTree::Branch()
1045 while ((obj = next())) {
1046 branchname.Form("%s", obj->GetName());
1047 if ((kSplitlevel > 1) && !obj->InheritsFrom(TClonesArray::Class())) {
1050 tree->Bronch(branchname, obj->ClassName(), fESDObjects->GetObjectRef(obj),
1051 kBufsize, kSplitlevel - 1);
1057 void AliESDEvent::ReadFromTree(TTree *tree){
1059 // Connect the ESDEvent to a tree
1062 Printf("%s %d AliESDEvent::ReadFromTree() Zero Pointer to Tree \n",(char*)__FILE__,__LINE__);
1066 if(!tree->GetTree())tree->LoadTree(0);
1068 // if we find the "ESD" branch on the tree we do have the old structure
1069 if(tree->GetBranch("ESD")) {
1070 char ** address = (char **)(tree->GetBranch("ESD")->GetAddress());
1071 // do we have the friend branch
1072 TBranch * esdFB = tree->GetBranch("ESDfriend.");
1073 char ** addressF = 0;
1074 if(esdFB)addressF = (char **)(esdFB->GetAddress());
1076 printf("%s %d AliESDEvent::ReadFromTree() Reading old Tree \n",(char*)__FILE__,__LINE__);
1077 tree->SetBranchAddress("ESD", &fESDOld);
1079 tree->SetBranchAddress("ESDfriend.",&fESDFriendOld);
1082 printf("%s %d AliESDEvent::ReadFromTree() Reading old Tree \n",(char*)__FILE__,__LINE__);
1083 printf("%s %d Branch already connected. Using existing branch address. \n",(char*)__FILE__,__LINE__);
1084 fESDOld = (AliESD*) (*address);
1085 // addressF can still be 0, since branch needs to switched on
1086 if(addressF)fESDFriendOld = (AliESDfriend*) (*addressF);
1089 // have already connected the old ESD structure... ?
1090 // reuse also the pointer of the AlliESDEvent
1091 // otherwise create new ones
1092 TList* connectedList = (TList*) (tree->GetUserInfo()->FindObject("ESDObjectsConnectedToTree"));
1095 // If connected use the connected list of objects
1096 if(fESDObjects!= connectedList){
1097 // protect when called twice
1098 fESDObjects->Delete();
1099 fESDObjects = connectedList;
1104 // The pointer to the friend changes when called twice via InitIO
1105 // since AliESDEvent is deleted
1106 TObject* oldf = FindListObject("AliESDfriend");
1109 newf = (TObject*)*addressF;
1111 if(newf!=0&&oldf!=newf){
1112 // remove the old reference
1113 // Should we also delete it? Or is this handled in TTree I/O
1114 // since it is created by the first SetBranchAddress
1115 fESDObjects->Remove(oldf);
1117 fESDObjects->Add(newf);
1124 CreateStdContent(); // create for copy
1125 // if we have the esdfriend add it, so we always can access it via the userinfo
1126 if(fESDFriendOld)AddObject(fESDFriendOld);
1127 // we are not owner of the list objects
1128 // must not delete it
1129 fESDObjects->SetOwner(kFALSE);
1130 fESDObjects->SetName("ESDObjectsConnectedToTree");
1131 tree->GetUserInfo()->Add(fESDObjects);
1138 // Try to find AliESDEvent
1139 AliESDEvent *esdEvent = 0;
1140 esdEvent = (AliESDEvent*)tree->GetTree()->GetUserInfo()->FindObject("AliESDEvent");
1142 // Check if already connected to tree
1143 TList* connectedList = (TList*) (tree->GetUserInfo()->FindObject("ESDObjectsConnectedToTree"));
1144 if (connectedList) {
1145 // If connected use the connected list if objects
1146 fESDObjects->Delete();
1147 fESDObjects = connectedList;
1153 // prevent a memory leak when reading back the TList
1156 // create a new TList from the UserInfo TList...
1157 // copy constructor does not work...
1158 fESDObjects = (TList*)(esdEvent->GetList()->Clone());
1159 fESDObjects->SetOwner(kFALSE);
1160 if(fESDObjects->GetEntries()<kESDListN){
1161 printf("%s %d AliESDEvent::ReadFromTree() TList contains less than the standard contents %d < %d \n",
1162 (char*)__FILE__,__LINE__,fESDObjects->GetEntries(),kESDListN);
1164 // set the branch addresses
1165 TIter next(fESDObjects);
1167 while((el=(TNamed*)next())){
1168 TString bname(el->GetName());
1169 if(bname.CompareTo("AliESDfriend")==0)
1171 // AliESDfriend does not have a name ...
1172 tree->SetBranchAddress("ESDfriend.",fESDObjects->GetObjectRef(el));
1175 // check if branch exists under this Name
1176 TBranch *br = tree->GetBranch(bname.Data());
1178 tree->SetBranchAddress(bname.Data(),fESDObjects->GetObjectRef(el));
1181 br = tree->GetBranch(Form("%s.",bname.Data()));
1183 tree->SetBranchAddress(Form("%s.",bname.Data()),fESDObjects->GetObjectRef(el));
1186 printf("%s %d AliESDEvent::ReadFromTree() No Branch found with Name %s or %s. \n",
1187 (char*)__FILE__,__LINE__,bname.Data(),bname.Data());
1194 // when reading back we are not owner of the list
1195 // must not delete it
1196 fESDObjects->SetOwner(kFALSE);
1197 fESDObjects->SetName("ESDObjectsConnectedToTree");
1198 // we are not owner of the list objects
1199 // must not delete it
1200 tree->GetUserInfo()->Add(fESDObjects);
1204 // we can't get the list from the user data, create standard content
1205 // and set it by hand (no ESDfriend at the moment
1207 TIter next(fESDObjects);
1209 while((el=(TNamed*)next())){
1210 TString bname(el->GetName());
1211 tree->SetBranchAddress(bname.Data(),fESDObjects->GetObjectRef(el));
1214 // when reading back we are not owner of the list
1215 // must not delete it
1216 fESDObjects->SetOwner(kFALSE);
1221 void AliESDEvent::CopyFromOldESD()
1223 // Method which copies over everthing from the old esd structure to the
1228 SetRunNumber(fESDOld->GetRunNumber());
1229 SetPeriodNumber(fESDOld->GetPeriodNumber());
1230 SetMagneticField(fESDOld->GetMagneticField());
1232 // leave out diamond ...
1233 // SetDiamond(const AliESDVertex *vertex) { fESDRun->SetDiamond(vertex);}
1236 SetTriggerMask(fESDOld->GetTriggerMask());
1237 SetOrbitNumber(fESDOld->GetOrbitNumber());
1238 SetTimeStamp(fESDOld->GetTimeStamp());
1239 SetEventType(fESDOld->GetEventType());
1240 SetEventNumberInFile(fESDOld->GetEventNumberInFile());
1241 SetBunchCrossNumber(fESDOld->GetBunchCrossNumber());
1242 SetTriggerCluster(fESDOld->GetTriggerCluster());
1246 SetZDC(fESDOld->GetZDCN1Energy(),
1247 fESDOld->GetZDCP1Energy(),
1248 fESDOld->GetZDCEMEnergy(),
1250 fESDOld->GetZDCN2Energy(),
1251 fESDOld->GetZDCP2Energy(),
1252 fESDOld->GetZDCParticipants());
1256 if(fESDOld->GetFMDData())SetFMDData(fESDOld->GetFMDData());
1260 SetT0zVertex(fESDOld->GetT0zVertex());
1261 SetT0(fESDOld->GetT0());
1265 if (fESDOld->GetVZEROData()) SetVZEROData(fESDOld->GetVZEROData());
1267 if(fESDOld->GetVertex())SetPrimaryVertexSPD(fESDOld->GetVertex());
1269 if(fESDOld->GetPrimaryVertex())SetPrimaryVertex(fESDOld->GetPrimaryVertex());
1271 if(fESDOld->GetMultiplicity())SetMultiplicity(fESDOld->GetMultiplicity());
1273 for(int i = 0;i<fESDOld->GetNumberOfTracks();i++){
1274 AddTrack(fESDOld->GetTrack(i));
1277 for(int i = 0;i<fESDOld->GetNumberOfMuonTracks();i++){
1278 AddMuonTrack(fESDOld->GetMuonTrack(i));
1281 for(int i = 0;i<fESDOld->GetNumberOfPmdTracks();i++){
1282 AddPmdTrack(fESDOld->GetPmdTrack(i));
1285 for(int i = 0;i<fESDOld->GetNumberOfTrdTracks();i++){
1286 AddTrdTrack(fESDOld->GetTrdTrack(i));
1289 for(int i = 0;i<fESDOld->GetNumberOfV0s();i++){
1290 AddV0(fESDOld->GetV0(i));
1293 for(int i = 0;i<fESDOld->GetNumberOfCascades();i++){
1294 AddCascade(fESDOld->GetCascade(i));
1297 for(int i = 0;i<fESDOld->GetNumberOfKinks();i++){
1298 AddKink(fESDOld->GetKink(i));
1302 for(int i = 0;i<fESDOld->GetNumberOfCaloClusters();i++){
1303 AddCaloCluster(fESDOld->GetCaloCluster(i));