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 "AliVertexerTracks.h"
53 #include "AliESDcascade.h"
54 #include "AliESDkink.h"
55 #include "AliESDtrack.h"
56 #include "AliESDHLTtrack.h"
57 #include "AliESDCaloCluster.h"
58 #include "AliESDCaloCells.h"
60 #include "AliESDFMD.h"
61 #include "AliESDVZERO.h"
62 #include "AliMultiplicity.h"
63 #include "AliRawDataErrorLog.h"
65 #include "AliESDACORDE.h"
70 // here we define the names, some classes are no TNamed, therefore the classnames
72 const char* AliESDEvent::fgkESDListName[kESDListN] = {"AliESDRun",
96 "AliRawDataErrorLogs",
99 //______________________________________________________________________________
100 AliESDEvent::AliESDEvent():
102 fESDObjects(new TList()),
116 fSPDPileupVertices(0),
117 fTrkPileupVertices(0),
126 fEMCALCells(0), fPHOSCells(0),
133 fFirstEMCALCluster(-1),
135 fFirstPHOSCluster(-1)
138 //______________________________________________________________________________
139 AliESDEvent::AliESDEvent(const AliESDEvent& esd):
141 fESDObjects(new TList()),
142 fESDRun(new AliESDRun(*esd.fESDRun)),
143 fHeader(new AliESDHeader(*esd.fHeader)),
144 fESDZDC(new AliESDZDC(*esd.fESDZDC)),
145 fESDFMD(new AliESDFMD(*esd.fESDFMD)),
146 fESDVZERO(new AliESDVZERO(*esd.fESDVZERO)),
147 fESDTZERO(new AliESDTZERO(*esd.fESDTZERO)),
148 fTPCVertex(new AliESDVertex(*esd.fTPCVertex)),
149 fSPDVertex(new AliESDVertex(*esd.fSPDVertex)),
150 fPrimaryVertex(new AliESDVertex(*esd.fPrimaryVertex)),
151 fSPDMult(new AliMultiplicity(*esd.fSPDMult)),
152 fPHOSTrigger(new AliESDCaloTrigger(*esd.fPHOSTrigger)),
153 fEMCALTrigger(new AliESDCaloTrigger(*esd.fEMCALTrigger)),
154 fESDACORDE(new AliESDACORDE(*esd.fESDACORDE)),
155 fSPDPileupVertices(new TClonesArray(*esd.fSPDPileupVertices)),
156 fTrkPileupVertices(new TClonesArray(*esd.fTrkPileupVertices)),
157 fTracks(new TClonesArray(*esd.fTracks)),
158 fMuonTracks(new TClonesArray(*esd.fMuonTracks)),
159 fPmdTracks(new TClonesArray(*esd.fPmdTracks)),
160 fTrdTracks(new TClonesArray(*esd.fTrdTracks)),
161 fV0s(new TClonesArray(*esd.fV0s)),
162 fCascades(new TClonesArray(*esd.fCascades)),
163 fKinks(new TClonesArray(*esd.fKinks)),
164 fCaloClusters(new TClonesArray(*esd.fCaloClusters)),
165 fEMCALCells(new AliESDCaloCells(*esd.fEMCALCells)),
166 fPHOSCells(new AliESDCaloCells(*esd.fPHOSCells)),
167 fErrorLogs(new TClonesArray(*esd.fErrorLogs)),
168 fESDOld(new AliESD(*esd.fESDOld)),
169 fESDFriendOld(new AliESDfriend(*esd.fESDFriendOld)),
170 fConnected(esd.fConnected),
171 fUseOwnList(esd.fUseOwnList),
172 fEMCALClusters(esd.fEMCALClusters),
173 fFirstEMCALCluster(esd.fFirstEMCALCluster),
174 fPHOSClusters(esd.fPHOSClusters),
175 fFirstPHOSCluster(esd.fFirstPHOSCluster)
178 // CKB init in the constructor list and only add here ...
183 AddObject(fESDVZERO);
184 AddObject(fESDTZERO);
185 AddObject(fTPCVertex);
186 AddObject(fSPDVertex);
187 AddObject(fPrimaryVertex);
189 AddObject(fPHOSTrigger);
190 AddObject(fEMCALTrigger);
191 AddObject(fSPDPileupVertices);
192 AddObject(fTrkPileupVertices);
194 AddObject(fMuonTracks);
195 AddObject(fPmdTracks);
196 AddObject(fTrdTracks);
198 AddObject(fCascades);
200 AddObject(fCaloClusters);
201 AddObject(fEMCALCells);
202 AddObject(fPHOSCells);
203 AddObject(fErrorLogs);
204 AddObject(fESDACORDE);
210 //______________________________________________________________________________
211 AliESDEvent & AliESDEvent::operator=(const AliESDEvent& source) {
213 // Assignment operator
215 if(&source == this) return *this;
216 AliVEvent::operator=(source);
218 // This assumes that the list is already created
219 // and that the virtual void Copy(Tobject&) function
220 // is correctly implemented in the derived class
221 // otherwise only TObject::Copy() will be used
225 if((fESDObjects->GetSize()==0)&&(source.fESDObjects->GetSize()>=kESDListN)){
226 // We cover the case that we do not yet have the
227 // standard content but the source has it
231 TIter next(source.GetList());
234 while ((its = next())) {
235 name.Form("%s", its->GetName());
236 TObject *mine = fESDObjects->FindObject(name.Data());
238 TClass* pClass=TClass::GetClass(its->ClassName());
240 AliWarning(Form("Can not find class description for entry %s (%s)\n",
241 its->ClassName(), name.Data()));
245 mine=(TObject*)pClass->New();
247 // not in this: can be added to list
248 AliWarning(Form("%s:%d Could not find %s for copying \n",
249 (char*)__FILE__,__LINE__,name.Data()));
252 if(mine->InheritsFrom("TNamed")){
253 ((TNamed*)mine)->SetName(name);
255 else if(mine->InheritsFrom("TCollection")){
256 if(mine->InheritsFrom("TClonesArray"))
257 dynamic_cast<TClonesArray*>(mine)->SetClass(dynamic_cast<TClonesArray*>(its)->GetClass());
258 dynamic_cast<TCollection*>(mine)->SetName(name);
260 AliDebug(1, Form("adding object %s of type %s", mine->GetName(), mine->ClassName()));
264 if(!its->InheritsFrom("TCollection")){
268 else if(its->InheritsFrom("TClonesArray")){
269 // Create or expand the tclonesarray pointers
270 // so we can directly copy to the object
271 TClonesArray *its_tca = (TClonesArray*)its;
272 TClonesArray *mine_tca = (TClonesArray*)mine;
274 // this leaves the capacity of the TClonesArray the same
275 // except for a factor of 2 increase when size > capacity
276 // does not release any memory occupied by the tca
277 mine_tca->ExpandCreate(its_tca->GetEntriesFast());
278 for(int i = 0;i < its_tca->GetEntriesFast();++i){
280 TObject *mine_tca_obj = mine_tca->At(i);
281 TObject *its_tca_obj = its_tca->At(i);
282 // no need to delete first
283 // pointers within the class should be handled by Copy()...
284 // Can there be Empty slots?
285 its_tca_obj->Copy(*mine_tca_obj);
289 AliWarning(Form("%s:%d cannot copy TCollection \n",
290 (char*)__FILE__,__LINE__));
294 fConnected = source.fConnected;
295 fUseOwnList = source.fUseOwnList;
296 fEMCALClusters = source.fEMCALClusters;
297 fFirstEMCALCluster = source.fFirstEMCALCluster;
298 fPHOSClusters = source.fPHOSClusters;
299 fFirstPHOSCluster = source.fFirstPHOSCluster;
307 //______________________________________________________________________________
308 AliESDEvent::~AliESDEvent()
311 // Standard destructor
314 // everthing on the list gets deleted automatically
317 if(fESDObjects&&!fConnected)
326 void AliESDEvent::Copy(TObject &obj) const {
328 // interface to TOBject::Copy
329 // Copies the content of this into obj!
330 // bascially obj = *this
332 if(this==&obj)return;
333 AliESDEvent *robj = dynamic_cast<AliESDEvent*>(&obj);
334 if(!robj)return; // not an AliESEvent
339 //______________________________________________________________________________
340 void AliESDEvent::Reset()
344 // Std content + Non std content
346 // Reset the standard contents
349 // reset for the old data without AliESDEvent...
350 if(fESDOld)fESDOld->Reset();
352 fESDFriendOld->~AliESDfriend();
353 new (fESDFriendOld) AliESDfriend();
357 if(fESDObjects->GetSize()>kESDListN){
358 // we have non std content
359 // this also covers esdfriends
360 for(int i = kESDListN;i < fESDObjects->GetSize();++i){
361 TObject *pObject = fESDObjects->At(i);
363 if(pObject->InheritsFrom(TClonesArray::Class())){
364 ((TClonesArray*)pObject)->Delete();
366 else if(!pObject->InheritsFrom(TCollection::Class())){
367 ResetWithPlacementNew(pObject);
370 AliWarning(Form("No reset for %s (%s)\n",
371 pObject->ClassName()));
378 Bool_t AliESDEvent::ResetWithPlacementNew(TObject *pObject){
379 Long_t dtoronly = TObject::GetDtorOnly();
380 TClass *pClass = TClass::GetClass(pObject->ClassName());
381 TObject::SetDtorOnly(pObject);
383 // Recreate with placement new
384 pClass->New(pObject);
385 // Restore the state.
386 TObject::SetDtorOnly((void*)dtoronly);
390 void AliESDEvent::ResetStdContent()
392 // Reset the standard contents
393 if(fESDRun) fESDRun->Reset();
394 if(fHeader) fHeader->Reset();
395 if(fESDZDC) fESDZDC->Reset();
400 // reset by callin d'to /c'tor keep the pointer
401 fESDVZERO->~AliESDVZERO();
402 new (fESDVZERO) AliESDVZERO();
405 fESDACORDE->~AliESDACORDE();
406 new (fESDACORDE) AliESDACORDE();
408 if(fESDTZERO) fESDTZERO->Reset();
409 // CKB no clear/reset implemented
411 fTPCVertex->~AliESDVertex();
412 new (fTPCVertex) AliESDVertex();
413 fTPCVertex->SetName(fgkESDListName[kTPCVertex]);
416 fSPDVertex->~AliESDVertex();
417 new (fSPDVertex) AliESDVertex();
418 fSPDVertex->SetName(fgkESDListName[kSPDVertex]);
421 fPrimaryVertex->~AliESDVertex();
422 new (fPrimaryVertex) AliESDVertex();
423 fPrimaryVertex->SetName(fgkESDListName[kPrimaryVertex]);
426 fSPDMult->~AliMultiplicity();
427 new (fSPDMult) AliMultiplicity();
429 if(fPHOSTrigger)fPHOSTrigger->Reset();
430 if(fEMCALTrigger)fEMCALTrigger->Reset();
431 if(fSPDPileupVertices)fSPDPileupVertices->Delete();
432 if(fTrkPileupVertices)fTrkPileupVertices->Delete();
433 if(fTracks)fTracks->Delete();
434 if(fMuonTracks)fMuonTracks->Delete();
435 if(fPmdTracks)fPmdTracks->Delete();
436 if(fTrdTracks)fTrdTracks->Delete();
437 if(fV0s)fV0s->Delete();
438 if(fCascades)fCascades->Delete();
439 if(fKinks)fKinks->Delete();
440 if(fCaloClusters)fCaloClusters->Delete();
441 if(fPHOSCells)fPHOSCells->DeleteContainer();
442 if(fEMCALCells)fEMCALCells->DeleteContainer();
443 if(fErrorLogs) fErrorLogs->Delete();
445 // don't reset fconnected fConnected and the list
448 fFirstEMCALCluster=-1;
450 fFirstPHOSCluster=-1;
454 Int_t AliESDEvent::AddV0(const AliESDv0 *v) {
458 TClonesArray &fv = *fV0s;
459 Int_t idx=fV0s->GetEntriesFast();
460 new(fv[idx]) AliESDv0(*v);
464 //______________________________________________________________________________
465 void AliESDEvent::Print(Option_t *) const
468 // Print header information of the event
470 printf("ESD run information\n");
471 printf("Event # in file %d Bunch crossing # %d Orbit # %d Period # %d Run # %d Trigger %lld Magnetic field %f \n",
472 GetEventNumberInFile(),
473 GetBunchCrossNumber(),
478 GetMagneticField() );
479 printf("Vertex: (%.4f +- %.4f, %.4f +- %.4f, %.4f +- %.4f) cm\n",
480 fPrimaryVertex->GetXv(), fPrimaryVertex->GetXRes(),
481 fPrimaryVertex->GetYv(), fPrimaryVertex->GetYRes(),
482 fPrimaryVertex->GetZv(), fPrimaryVertex->GetZRes());
483 printf("Mean vertex in RUN: X=%.4f Y=%.4f cm\n",
484 GetDiamondX(),GetDiamondY());
485 printf("SPD Multiplicity. Number of tracklets %d \n",
486 fSPDMult->GetNumberOfTracklets());
487 printf("Number of pileup primary vertices reconstructed with SPD %d\n",
488 GetNumberOfPileupVerticesSPD());
489 printf("Number of pileup primary vertices reconstructed using the tracks %d\n",
490 GetNumberOfPileupVerticesTracks());
491 printf("Number of tracks: \n");
492 printf(" charged %d\n", GetNumberOfTracks());
493 printf(" muon %d\n", GetNumberOfMuonTracks());
494 printf(" pmd %d\n", GetNumberOfPmdTracks());
495 printf(" trd %d\n", GetNumberOfTrdTracks());
496 printf(" v0 %d\n", GetNumberOfV0s());
497 printf(" cascades %d\n", GetNumberOfCascades());
498 printf(" kinks %d\n", GetNumberOfKinks());
499 if(fPHOSCells)printf(" PHOSCells %d\n", fPHOSCells->GetNumberOfCells());
500 else printf(" PHOSCells not in the Event\n");
501 if(fEMCALCells)printf(" EMCALCells %d\n", fEMCALCells->GetNumberOfCells());
502 else printf(" EMCALCells not in the Event\n");
503 printf(" CaloClusters %d\n", GetNumberOfCaloClusters());
504 printf(" phos %d\n", GetNumberOfPHOSClusters());
505 printf(" emcal %d\n", GetNumberOfEMCALClusters());
506 printf(" FMD %s\n", (fESDFMD ? "yes" : "no"));
507 printf(" VZERO %s\n", (fESDVZERO ? "yes" : "no"));
512 void AliESDEvent::SetESDfriend(const AliESDfriend *ev) const {
514 // Attaches the complementary info to the ESD
518 // to be sure that we set the tracks also
519 // in case of old esds
520 // if(fESDOld)CopyFromOldESD();
522 Int_t ntrk=ev->GetNumberOfTracks();
524 for (Int_t i=0; i<ntrk; i++) {
525 const AliESDfriendTrack *f=ev->GetTrack(i);
526 GetTrack(i)->SetFriendTrack(f);
530 Bool_t AliESDEvent::RemoveKink(Int_t rm) const {
531 // ---------------------------------------------------------
532 // Remove a kink candidate and references to it from ESD,
533 // if this candidate does not come from a reconstructed decay
534 // Not yet implemented...
535 // ---------------------------------------------------------
536 Int_t last=GetNumberOfKinks()-1;
537 if ((rm<0)||(rm>last)) return kFALSE;
542 Bool_t AliESDEvent::RemoveV0(Int_t rm) const {
543 // ---------------------------------------------------------
544 // Remove a V0 candidate and references to it from ESD,
545 // if this candidate does not come from a reconstructed decay
546 // ---------------------------------------------------------
547 Int_t last=GetNumberOfV0s()-1;
548 if ((rm<0)||(rm>last)) return kFALSE;
550 AliESDv0 *v0=GetV0(rm);
551 Int_t idxP=v0->GetPindex(), idxN=v0->GetNindex();
554 Int_t lastIdxP=v0->GetPindex(), lastIdxN=v0->GetNindex();
558 // Check if this V0 comes from a reconstructed decay
559 Int_t ncs=GetNumberOfCascades();
560 for (Int_t n=0; n<ncs; n++) {
561 AliESDcascade *cs=GetCascade(n);
563 Int_t csIdxP=cs->GetPindex();
564 Int_t csIdxN=cs->GetNindex();
567 if (idxN==csIdxN) return kFALSE;
569 if (csIdxP==lastIdxP)
570 if (csIdxN==lastIdxN) used++;
573 //Replace the removed V0 with the last V0
574 TClonesArray &a=*fV0s;
575 delete a.RemoveAt(rm);
577 if (rm==last) return kTRUE;
579 //v0 is pointing to the last V0 candidate...
580 new (a[rm]) AliESDv0(*v0);
581 delete a.RemoveAt(last);
583 if (!used) return kTRUE;
586 // Remap the indices of the daughters of reconstructed decays
587 for (Int_t n=0; n<ncs; n++) {
588 AliESDcascade *cs=GetCascade(n);
591 Int_t csIdxP=cs->GetPindex();
592 Int_t csIdxN=cs->GetNindex();
594 if (csIdxP==lastIdxP)
595 if (csIdxN==lastIdxN) {
596 cs->AliESDv0::SetIndex(1,idxP);
597 cs->AliESDv0::SetIndex(0,idxN);
599 if (!used) return kTRUE;
606 Bool_t AliESDEvent::RemoveTrack(Int_t rm) const {
607 // ---------------------------------------------------------
608 // Remove a track and references to it from ESD,
609 // if this track does not come from a reconstructed decay
610 // ---------------------------------------------------------
611 Int_t last=GetNumberOfTracks()-1;
612 if ((rm<0)||(rm>last)) return kFALSE;
616 // Check if this track comes from the reconstructed primary vertices
617 if (fTPCVertex && fTPCVertex->GetStatus()) {
618 UShort_t *primIdx=fTPCVertex->GetIndices();
619 Int_t n=fTPCVertex->GetNIndices();
621 Int_t idx=Int_t(primIdx[n]);
622 if (rm==idx) return kFALSE;
623 if (idx==last) used++;
626 if (fPrimaryVertex && fPrimaryVertex->GetStatus()) {
627 UShort_t *primIdx=fPrimaryVertex->GetIndices();
628 Int_t n=fPrimaryVertex->GetNIndices();
630 Int_t idx=Int_t(primIdx[n]);
631 if (rm==idx) return kFALSE;
632 if (idx==last) used++;
636 // Check if this track comes from a reconstructed decay
637 Int_t nv0=GetNumberOfV0s();
638 for (Int_t n=0; n<nv0; n++) {
639 AliESDv0 *v0=GetV0(n);
641 Int_t idx=v0->GetNindex();
642 if (rm==idx) return kFALSE;
643 if (idx==last) used++;
646 if (rm==idx) return kFALSE;
647 if (idx==last) used++;
650 Int_t ncs=GetNumberOfCascades();
651 for (Int_t n=0; n<ncs; n++) {
652 AliESDcascade *cs=GetCascade(n);
654 Int_t idx=cs->GetIndex();
655 if (rm==idx) return kFALSE;
656 if (idx==last) used++;
660 if (rm==idx) return kFALSE;
661 if (idx==last) used++;
664 if (rm==idx) return kFALSE;
665 if (idx==last) used++;
668 Int_t nkn=GetNumberOfKinks();
669 for (Int_t n=0; n<nkn; n++) {
670 AliESDkink *kn=GetKink(n);
672 Int_t idx=kn->GetIndex(0);
673 if (rm==idx) return kFALSE;
674 if (idx==last) used++;
677 if (rm==idx) return kFALSE;
678 if (idx==last) used++;
681 // Check if this track is associated with a CaloCluster
682 Int_t ncl=GetNumberOfCaloClusters();
683 for (Int_t n=0; n<ncl; n++) {
684 AliESDCaloCluster *cluster=GetCaloCluster(n);
685 TArrayI *arr=cluster->GetTracksMatched();
686 Int_t s=arr->GetSize();
688 Int_t idx=arr->At(s);
689 if (rm==idx) return kFALSE;
690 if (idx==last) used++;
696 //Replace the removed track with the last track
697 TClonesArray &a=*fTracks;
698 delete a.RemoveAt(rm);
700 if (rm==last) return kTRUE;
702 AliESDtrack *t=GetTrack(last);
704 new (a[rm]) AliESDtrack(*t);
705 delete a.RemoveAt(last);
708 if (!used) return kTRUE;
711 // Remap the indices of the tracks used for the primary vertex reconstruction
712 if (fTPCVertex && fTPCVertex->GetStatus()) {
713 UShort_t *primIdx=fTPCVertex->GetIndices();
714 Int_t n=fTPCVertex->GetNIndices();
716 Int_t idx=Int_t(primIdx[n]);
718 primIdx[n]=Short_t(rm);
720 if (!used) return kTRUE;
724 if (fPrimaryVertex && fPrimaryVertex->GetStatus()) {
725 UShort_t *primIdx=fPrimaryVertex->GetIndices();
726 Int_t n=fPrimaryVertex->GetNIndices();
728 Int_t idx=Int_t(primIdx[n]);
730 primIdx[n]=Short_t(rm);
732 if (!used) return kTRUE;
737 // Remap the indices of the daughters of reconstructed decays
738 for (Int_t n=0; n<nv0; n++) {
739 AliESDv0 *v0=GetV0(n);
740 if (v0->GetIndex(0)==last) {
743 if (!used) return kTRUE;
745 if (v0->GetIndex(1)==last) {
748 if (!used) return kTRUE;
752 for (Int_t n=0; n<ncs; n++) {
753 AliESDcascade *cs=GetCascade(n);
754 if (cs->GetIndex()==last) {
757 if (!used) return kTRUE;
760 if (v0->GetIndex(0)==last) {
763 if (!used) return kTRUE;
765 if (v0->GetIndex(1)==last) {
768 if (!used) return kTRUE;
772 for (Int_t n=0; n<nkn; n++) {
773 AliESDkink *kn=GetKink(n);
774 if (kn->GetIndex(0)==last) {
777 if (!used) return kTRUE;
779 if (kn->GetIndex(1)==last) {
782 if (!used) return kTRUE;
786 // Remap the indices of the tracks accosicated with CaloClusters
787 for (Int_t n=0; n<ncl; n++) {
788 AliESDCaloCluster *cluster=GetCaloCluster(n);
789 TArrayI *arr=cluster->GetTracksMatched();
790 Int_t s=arr->GetSize();
792 Int_t idx=arr->At(s);
796 if (!used) return kTRUE;
805 Bool_t AliESDEvent::Clean(Float_t *cleanPars) {
807 // Remove the data which are not needed for the physics analysis.
809 // 1) Cleaning the V0 candidates
810 // ---------------------------
811 // If the cosine of the V0 pointing angle "csp" and
812 // the DCA between the daughter tracks "dca" does not satisfy
815 // csp > cleanPars[1] + dca/cleanPars[0]*(1.- cleanPars[1])
817 // an attempt to remove this V0 candidate from ESD is made.
819 // The V0 candidate gets removed if it does not belong to any
820 // recosntructed cascade decay
822 // 12.11.2007, optimal values: cleanPars[0]=0.5, cleanPars[1]=0.999
824 // 2) Cleaning the tracks
825 // ----------------------
826 // If track's transverse parameter is larger than cleanPars[2]
828 // track's longitudinal parameter is larger than cleanPars[3]
829 // an attempt to remove this track from ESD is made.
831 // The track gets removed if it does not come
832 // from a reconstructed decay
836 Float_t dcaMax=cleanPars[0];
837 Float_t cspMin=cleanPars[1];
839 Int_t nV0s=GetNumberOfV0s();
840 for (Int_t i=nV0s-1; i>=0; i--) {
841 AliESDv0 *v0=GetV0(i);
843 Float_t dca=v0->GetDcaV0Daughters();
844 Float_t csp=v0->GetV0CosineOfPointingAngle();
845 Float_t cspcut=cspMin + dca/dcaMax*(1.-cspMin);
846 if (csp > cspcut) continue;
847 if (RemoveV0(i)) rc=kTRUE;
851 Float_t dmax=cleanPars[2], zmax=cleanPars[3];
853 const AliESDVertex *vertex=GetPrimaryVertexSPD();
854 Bool_t vtxOK=vertex->GetStatus();
856 Int_t nTracks=GetNumberOfTracks();
857 for (Int_t i=nTracks-1; i>=0; i--) {
858 AliESDtrack *track=GetTrack(i);
859 Float_t xy,z; track->GetImpactParameters(xy,z);
860 if ((TMath::Abs(xy) > dmax) || (vtxOK && (TMath::Abs(z) > zmax))) {
861 if (RemoveTrack(i)) rc=kTRUE;
868 Char_t AliESDEvent::AddPileupVertexSPD(const AliESDVertex *vtx)
870 // Add a pileup primary vertex reconstructed with SPD
871 TClonesArray &ftr = *fSPDPileupVertices;
872 Char_t n=Char_t(ftr.GetEntriesFast());
873 AliESDVertex *vertex = new(ftr[n]) AliESDVertex(*vtx);
878 Char_t AliESDEvent::AddPileupVertexTracks(const AliESDVertex *vtx)
880 // Add a pileup primary vertex reconstructed with SPD
881 TClonesArray &ftr = *fTrkPileupVertices;
882 Char_t n=Char_t(ftr.GetEntriesFast());
883 AliESDVertex *vertex = new(ftr[n]) AliESDVertex(*vtx);
888 Int_t AliESDEvent::AddTrack(const AliESDtrack *t)
891 TClonesArray &ftr = *fTracks;
892 AliESDtrack * track = new(ftr[fTracks->GetEntriesFast()])AliESDtrack(*t);
893 track->SetID(fTracks->GetEntriesFast()-1);
894 return track->GetID();
897 void AliESDEvent::AddMuonTrack(const AliESDMuonTrack *t)
899 TClonesArray &fmu = *fMuonTracks;
900 new(fmu[fMuonTracks->GetEntriesFast()]) AliESDMuonTrack(*t);
903 void AliESDEvent::AddPmdTrack(const AliESDPmdTrack *t)
905 TClonesArray &fpmd = *fPmdTracks;
906 new(fpmd[fPmdTracks->GetEntriesFast()]) AliESDPmdTrack(*t);
909 void AliESDEvent::AddTrdTrack(const AliESDTrdTrack *t)
911 TClonesArray &ftrd = *fTrdTracks;
912 new(ftrd[fTrdTracks->GetEntriesFast()]) AliESDTrdTrack(*t);
918 Int_t AliESDEvent::AddKink(const AliESDkink *c)
921 TClonesArray &fk = *fKinks;
922 AliESDkink * kink = new(fk[fKinks->GetEntriesFast()]) AliESDkink(*c);
923 kink->SetID(fKinks->GetEntriesFast()); // CKB different from the other imps..
924 return fKinks->GetEntriesFast()-1;
928 void AliESDEvent::AddCascade(const AliESDcascade *c)
930 TClonesArray &fc = *fCascades;
931 new(fc[fCascades->GetEntriesFast()]) AliESDcascade(*c);
935 Int_t AliESDEvent::AddCaloCluster(const AliESDCaloCluster *c)
938 TClonesArray &fc = *fCaloClusters;
939 AliESDCaloCluster *clus = new(fc[fCaloClusters->GetEntriesFast()]) AliESDCaloCluster(*c);
940 clus->SetID(fCaloClusters->GetEntriesFast()-1);
941 return fCaloClusters->GetEntriesFast()-1;
945 void AliESDEvent::AddRawDataErrorLog(const AliRawDataErrorLog *log) const {
946 TClonesArray &errlogs = *fErrorLogs;
947 new(errlogs[errlogs.GetEntriesFast()]) AliRawDataErrorLog(*log);
950 void AliESDEvent::SetPrimaryVertexTPC(const AliESDVertex *vertex)
952 // Set the TPC vertex
953 // use already allocated space
955 *fTPCVertex = *vertex;
956 fTPCVertex->SetName(fgkESDListName[kTPCVertex]);
960 void AliESDEvent::SetPrimaryVertexSPD(const AliESDVertex *vertex)
962 // Set the SPD vertex
963 // use already allocated space
965 *fSPDVertex = *vertex;
966 fSPDVertex->SetName(fgkESDListName[kSPDVertex]);
970 void AliESDEvent::SetPrimaryVertexTracks(const AliESDVertex *vertex)
972 // Set the primary vertex reconstructed using he ESD tracks.
973 // use already allocated space
975 *fPrimaryVertex = *vertex;
976 fPrimaryVertex->SetName(fgkESDListName[kPrimaryVertex]);
980 const AliESDVertex * AliESDEvent::GetPrimaryVertex() const
983 // Get the "best" available reconstructed primary vertex.
986 if (fPrimaryVertex->GetStatus()) return fPrimaryVertex;
989 if (fSPDVertex->GetStatus()) return fSPDVertex;
991 if(fTPCVertex) return fTPCVertex;
993 AliWarning("No primary vertex available. Returning the \"default\"...");
997 AliESDVertex * AliESDEvent::PrimaryVertexTracksUnconstrained() const
1000 // Removes diamond constraint from fPrimaryVertex (reconstructed with tracks)
1001 // Returns a AliESDVertex which has to be deleted by the user
1003 if(!fPrimaryVertex) {
1004 AliWarning("No primary vertex from tracks available.");
1007 if(!fPrimaryVertex->GetStatus()) {
1008 AliWarning("No primary vertex from tracks available.");
1012 AliVertexerTracks vertexer(GetMagneticField());
1013 Float_t diamondxyz[3]={(Float_t)GetDiamondX(),(Float_t)GetDiamondY(),0.};
1014 Float_t diamondcovxy[3]; GetDiamondCovXY(diamondcovxy);
1015 Float_t diamondcov[6]={diamondcovxy[0],diamondcovxy[1],diamondcovxy[2],0.,0.,7.};
1016 AliESDVertex *vertex =
1017 (AliESDVertex*)vertexer.RemoveConstraintFromVertex(fPrimaryVertex,diamondxyz,diamondcov);
1022 void AliESDEvent::SetMultiplicity(const AliMultiplicity *mul)
1024 // Set the SPD Multiplicity
1031 void AliESDEvent::SetFMDData(AliESDFMD * obj)
1033 // use already allocated space
1039 void AliESDEvent::SetVZEROData(AliESDVZERO * obj)
1041 // use already allocated space
1046 void AliESDEvent::SetACORDEData(AliESDACORDE * obj)
1053 void AliESDEvent::GetESDfriend(AliESDfriend *ev) const
1056 // Extracts the complementary info from the ESD
1060 Int_t ntrk=GetNumberOfTracks();
1062 for (Int_t i=0; i<ntrk; i++) {
1063 AliESDtrack *t=GetTrack(i);
1064 const AliESDfriendTrack *f=t->GetFriendTrack();
1067 t->ReleaseESDfriendTrack();// Not to have two copies of "friendTrack"
1071 AliESDfriend *fr = (AliESDfriend*)(const_cast<AliESDEvent*>(this)->FindListObject("AliESDfriend"));
1072 if (fr) ev->SetVZEROfriend(fr->GetVZEROfriend());
1075 void AliESDEvent::AddObject(TObject* obj)
1077 // Add an object to the list of object.
1078 // Please be aware that in order to increase performance you should
1079 // refrain from using TObjArrays (if possible). Use TClonesArrays, instead.
1080 fESDObjects->SetOwner(kTRUE);
1081 fESDObjects->AddLast(obj);
1085 void AliESDEvent::GetStdContent()
1087 // set pointers for standard content
1088 // get by name much safer and not a big overhead since not called very often
1090 fESDRun = (AliESDRun*)fESDObjects->FindObject(fgkESDListName[kESDRun]);
1091 fHeader = (AliESDHeader*)fESDObjects->FindObject(fgkESDListName[kHeader]);
1092 fESDZDC = (AliESDZDC*)fESDObjects->FindObject(fgkESDListName[kESDZDC]);
1093 fESDFMD = (AliESDFMD*)fESDObjects->FindObject(fgkESDListName[kESDFMD]);
1094 fESDVZERO = (AliESDVZERO*)fESDObjects->FindObject(fgkESDListName[kESDVZERO]);
1095 fESDTZERO = (AliESDTZERO*)fESDObjects->FindObject(fgkESDListName[kESDTZERO]);
1096 fTPCVertex = (AliESDVertex*)fESDObjects->FindObject(fgkESDListName[kTPCVertex]);
1097 fSPDVertex = (AliESDVertex*)fESDObjects->FindObject(fgkESDListName[kSPDVertex]);
1098 fPrimaryVertex = (AliESDVertex*)fESDObjects->FindObject(fgkESDListName[kPrimaryVertex]);
1099 fSPDMult = (AliMultiplicity*)fESDObjects->FindObject(fgkESDListName[kSPDMult]);
1100 fPHOSTrigger = (AliESDCaloTrigger*)fESDObjects->FindObject(fgkESDListName[kPHOSTrigger]);
1101 fEMCALTrigger = (AliESDCaloTrigger*)fESDObjects->FindObject(fgkESDListName[kEMCALTrigger]);
1102 fSPDPileupVertices = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kSPDPileupVertices]);
1103 fTrkPileupVertices = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kTrkPileupVertices]);
1104 fTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kTracks]);
1105 fMuonTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kMuonTracks]);
1106 fPmdTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kPmdTracks]);
1107 fTrdTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kTrdTracks]);
1108 fV0s = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kV0s]);
1109 fCascades = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kCascades]);
1110 fKinks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kKinks]);
1111 fCaloClusters = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kCaloClusters]);
1112 fEMCALCells = (AliESDCaloCells*)fESDObjects->FindObject(fgkESDListName[kEMCALCells]);
1113 fPHOSCells = (AliESDCaloCells*)fESDObjects->FindObject(fgkESDListName[kPHOSCells]);
1114 fErrorLogs = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kErrorLogs]);
1115 fESDACORDE = (AliESDACORDE*)fESDObjects->FindObject(fgkESDListName[kESDACORDE]);
1119 void AliESDEvent::SetStdNames(){
1120 // Set the names of the standard contents
1122 if(fESDObjects->GetEntries()>=kESDListN){
1123 for(int i = 0;i < fESDObjects->GetEntries() && i<kESDListN;i++){
1124 TObject *fObj = fESDObjects->At(i);
1125 if(fObj->InheritsFrom("TNamed")){
1126 ((TNamed*)fObj)->SetName(fgkESDListName[i]);
1128 else if(fObj->InheritsFrom("TClonesArray")){
1129 ((TClonesArray*)fObj)->SetName(fgkESDListName[i]);
1134 AliWarning("Std Entries missing");
1139 void AliESDEvent::CreateStdContent(Bool_t bUseThisList){
1140 fUseOwnList = bUseThisList;
1144 void AliESDEvent::CreateStdContent()
1146 // create the standard AOD content and set pointers
1148 // create standard objects and add them to the TList of objects
1149 AddObject(new AliESDRun());
1150 AddObject(new AliESDHeader());
1151 AddObject(new AliESDZDC());
1152 AddObject(new AliESDFMD());
1153 AddObject(new AliESDVZERO());
1154 AddObject(new AliESDTZERO());
1155 AddObject(new AliESDVertex());
1156 AddObject(new AliESDVertex());
1157 AddObject(new AliESDVertex());
1158 AddObject(new AliMultiplicity());
1159 AddObject(new AliESDCaloTrigger());
1160 AddObject(new AliESDCaloTrigger());
1161 AddObject(new TClonesArray("AliESDVertex",0));
1162 AddObject(new TClonesArray("AliESDVertex",0));
1163 AddObject(new TClonesArray("AliESDtrack",0));
1164 AddObject(new TClonesArray("AliESDMuonTrack",0));
1165 AddObject(new TClonesArray("AliESDPmdTrack",0));
1166 AddObject(new TClonesArray("AliESDTrdTrack",0));
1167 AddObject(new TClonesArray("AliESDv0",0));
1168 AddObject(new TClonesArray("AliESDcascade",0));
1169 AddObject(new TClonesArray("AliESDkink",0));
1170 AddObject(new TClonesArray("AliESDCaloCluster",0));
1171 AddObject(new AliESDCaloCells());
1172 AddObject(new AliESDCaloCells());
1173 AddObject(new TClonesArray("AliRawDataErrorLog",0));
1174 AddObject(new AliESDACORDE());
1176 // check the order of the indices against enum...
1180 // read back pointers
1184 TObject* AliESDEvent::FindListObject(const char *name){
1186 // Find object with name "name" in the list of branches
1189 return fESDObjects->FindObject(name);
1194 Int_t AliESDEvent::GetPHOSClusters(TRefArray *clusters) const
1196 // fills the provided TRefArray with all found phos clusters
1200 AliESDCaloCluster *cl = 0;
1201 for (Int_t i = 0; i < GetNumberOfCaloClusters(); i++) {
1203 if ( (cl = GetCaloCluster(i)) ) {
1206 AliDebug(1,Form("IsPHOS cluster %d Size: %d \n",i,clusters->GetEntriesFast()));
1210 return clusters->GetEntriesFast();
1213 Int_t AliESDEvent::GetEMCALClusters(TRefArray *clusters) const
1215 // fills the provided TRefArray with all found emcal clusters
1219 AliESDCaloCluster *cl = 0;
1220 for (Int_t i = 0; i < GetNumberOfCaloClusters(); i++) {
1222 if ( (cl = GetCaloCluster(i)) ) {
1225 AliDebug(1,Form("IsEMCAL cluster %d Size: %d \n",i,clusters->GetEntriesFast()));
1229 return clusters->GetEntriesFast();
1232 void AliESDEvent::WriteToTree(TTree* tree) const {
1233 // Book the branches as in TTree::Branch(TCollection*)
1234 // but add a "." at the end of top level branches which are
1235 // not a TClonesArray
1239 TIter next(fESDObjects);
1240 const Int_t kSplitlevel = 99; // default value in TTree::Branch()
1241 const Int_t kBufsize = 32000; // default value in TTree::Branch()
1244 while ((obj = next())) {
1245 branchname.Form("%s", obj->GetName());
1246 if(branchname.CompareTo("AliESDfriend")==0)branchname = "ESDfriend.";
1247 if ((kSplitlevel > 1) && !obj->InheritsFrom(TClonesArray::Class())) {
1248 if(!branchname.EndsWith("."))branchname += ".";
1250 if (!tree->FindBranch(branchname)) {
1251 tree->Bronch(branchname, obj->ClassName(), fESDObjects->GetObjectRef(obj),
1252 kBufsize, kSplitlevel - 1);
1258 void AliESDEvent::ReadFromTree(TTree *tree, Option_t* opt){
1260 // Connect the ESDEvent to a tree
1263 AliWarning("AliESDEvent::ReadFromTree() Zero Pointer to Tree \n");
1267 if(!tree->GetTree())tree->LoadTree(0);
1269 // if we find the "ESD" branch on the tree we do have the old structure
1270 if(tree->GetBranch("ESD")) {
1271 char ** address = (char **)(tree->GetBranch("ESD")->GetAddress());
1272 // do we have the friend branch
1273 TBranch * esdFB = tree->GetBranch("ESDfriend.");
1274 char ** addressF = 0;
1275 if(esdFB)addressF = (char **)(esdFB->GetAddress());
1277 AliInfo("AliESDEvent::ReadFromTree() Reading old Tree");
1278 tree->SetBranchAddress("ESD", &fESDOld);
1280 tree->SetBranchAddress("ESDfriend.",&fESDFriendOld);
1283 AliInfo("AliESDEvent::ReadFromTree() Reading old Tree");
1284 AliInfo("Branch already connected. Using existing branch address.");
1285 fESDOld = (AliESD*) (*address);
1286 // addressF can still be 0, since branch needs to switched on
1287 if(addressF)fESDFriendOld = (AliESDfriend*) (*addressF);
1290 // have already connected the old ESD structure... ?
1291 // reuse also the pointer of the AlliESDEvent
1292 // otherwise create new ones
1293 TList* connectedList = (TList*) (tree->GetUserInfo()->FindObject("ESDObjectsConnectedToTree"));
1296 // If connected use the connected list of objects
1297 if(fESDObjects!= connectedList){
1298 // protect when called twice
1299 fESDObjects->Delete();
1300 fESDObjects = connectedList;
1305 // The pointer to the friend changes when called twice via InitIO
1306 // since AliESDEvent is deleted
1307 TObject* oldf = FindListObject("AliESDfriend");
1310 newf = (TObject*)*addressF;
1312 if(newf!=0&&oldf!=newf){
1313 // remove the old reference
1314 // Should we also delete it? Or is this handled in TTree I/O
1315 // since it is created by the first SetBranchAddress
1316 fESDObjects->Remove(oldf);
1318 fESDObjects->Add(newf);
1325 CreateStdContent(); // create for copy
1326 // if we have the esdfriend add it, so we always can access it via the userinfo
1327 if(fESDFriendOld)AddObject(fESDFriendOld);
1328 // we are not owner of the list objects
1329 // must not delete it
1330 fESDObjects->SetOwner(kTRUE);
1331 fESDObjects->SetName("ESDObjectsConnectedToTree");
1332 tree->GetUserInfo()->Add(fESDObjects);
1340 // Try to find AliESDEvent
1341 AliESDEvent *esdEvent = 0;
1342 esdEvent = (AliESDEvent*)tree->GetTree()->GetUserInfo()->FindObject("AliESDEvent");
1344 // Check if already connected to tree
1346 TList* connectedList = (TList*) (tree->GetUserInfo()->FindObject("ESDObjectsConnectedToTree"));
1349 if (connectedList && (strcmp(opt, "reconnect"))) {
1350 // If connected use the connected list if objects
1351 fESDObjects->Delete();
1352 fESDObjects = connectedList;
1359 // prevent a memory leak when reading back the TList
1360 if (!(strcmp(opt, "reconnect"))) fESDObjects->Delete();
1363 // create a new TList from the UserInfo TList...
1364 // copy constructor does not work...
1365 fESDObjects = (TList*)(esdEvent->GetList()->Clone());
1366 fESDObjects->SetOwner(kTRUE);
1368 else if ( fESDObjects->GetEntries()==0){
1369 // at least create the std content if we want to read to our list
1374 // we only need new things in the list if we do no already have it..
1375 // TODO just add new entries
1377 if(fESDObjects->GetEntries()<kESDListN){
1378 AliWarning(Form("AliESDEvent::ReadFromTree() TList contains less than the standard contents %d < %d \n",
1379 fESDObjects->GetEntries(),kESDListN));
1381 // set the branch addresses
1382 TIter next(fESDObjects);
1384 while((el=(TNamed*)next())){
1385 TString bname(el->GetName());
1386 if(bname.CompareTo("AliESDfriend")==0)
1388 // AliESDfriend does not have a name ...
1389 tree->SetBranchAddress("ESDfriend.",fESDObjects->GetObjectRef(el));
1392 // check if branch exists under this Name
1393 TBranch *br = tree->GetBranch(bname.Data());
1395 tree->SetBranchAddress(bname.Data(),fESDObjects->GetObjectRef(el));
1398 br = tree->GetBranch(Form("%s.",bname.Data()));
1400 tree->SetBranchAddress(Form("%s.",bname.Data()),fESDObjects->GetObjectRef(el));
1403 AliWarning(Form("AliESDEvent::ReadFromTree() No Branch found with Name %s or %s.",bname.Data(),bname.Data()));
1410 // when reading back we are not owner of the list
1411 // must not delete it
1412 fESDObjects->SetOwner(kTRUE);
1413 fESDObjects->SetName("ESDObjectsConnectedToTree");
1414 // we are not owner of the list objects
1415 // must not delete it
1416 tree->GetUserInfo()->Add(fESDObjects);
1417 tree->GetUserInfo()->SetOwner(kFALSE);
1421 // we can't get the list from the user data, create standard content
1422 // and set it by hand (no ESDfriend at the moment
1424 TIter next(fESDObjects);
1426 while((el=(TNamed*)next())){
1427 TString bname(el->GetName());
1428 TBranch *br = tree->GetBranch(bname.Data());
1430 tree->SetBranchAddress(bname.Data(),fESDObjects->GetObjectRef(el));
1433 br = tree->GetBranch(Form("%s.",bname.Data()));
1435 tree->SetBranchAddress(Form("%s.",bname.Data()),fESDObjects->GetObjectRef(el));
1440 // when reading back we are not owner of the list
1441 // must not delete it
1442 fESDObjects->SetOwner(kFALSE);
1447 void AliESDEvent::CopyFromOldESD()
1449 // Method which copies over everthing from the old esd structure to the
1454 SetRunNumber(fESDOld->GetRunNumber());
1455 SetPeriodNumber(fESDOld->GetPeriodNumber());
1456 SetMagneticField(fESDOld->GetMagneticField());
1458 // leave out diamond ...
1459 // SetDiamond(const AliESDVertex *vertex) { fESDRun->SetDiamond(vertex);}
1462 SetTriggerMask(fESDOld->GetTriggerMask());
1463 SetOrbitNumber(fESDOld->GetOrbitNumber());
1464 SetTimeStamp(fESDOld->GetTimeStamp());
1465 SetEventType(fESDOld->GetEventType());
1466 SetEventNumberInFile(fESDOld->GetEventNumberInFile());
1467 SetBunchCrossNumber(fESDOld->GetBunchCrossNumber());
1468 SetTriggerCluster(fESDOld->GetTriggerCluster());
1472 SetZDC(fESDOld->GetZDCN1Energy(),
1473 fESDOld->GetZDCP1Energy(),
1474 fESDOld->GetZDCEMEnergy(),
1476 fESDOld->GetZDCN2Energy(),
1477 fESDOld->GetZDCP2Energy(),
1478 fESDOld->GetZDCParticipants(),
1488 if(fESDOld->GetFMDData())SetFMDData(fESDOld->GetFMDData());
1492 SetT0zVertex(fESDOld->GetT0zVertex());
1493 SetT0(fESDOld->GetT0());
1497 if (fESDOld->GetVZEROData()) SetVZEROData(fESDOld->GetVZEROData());
1499 if(fESDOld->GetVertex())SetPrimaryVertexSPD(fESDOld->GetVertex());
1501 if(fESDOld->GetPrimaryVertex())SetPrimaryVertexTracks(fESDOld->GetPrimaryVertex());
1503 if(fESDOld->GetMultiplicity())SetMultiplicity(fESDOld->GetMultiplicity());
1505 for(int i = 0;i<fESDOld->GetNumberOfTracks();i++){
1506 AddTrack(fESDOld->GetTrack(i));
1509 for(int i = 0;i<fESDOld->GetNumberOfMuonTracks();i++){
1510 AddMuonTrack(fESDOld->GetMuonTrack(i));
1513 for(int i = 0;i<fESDOld->GetNumberOfPmdTracks();i++){
1514 AddPmdTrack(fESDOld->GetPmdTrack(i));
1517 for(int i = 0;i<fESDOld->GetNumberOfTrdTracks();i++){
1518 AddTrdTrack(fESDOld->GetTrdTrack(i));
1521 for(int i = 0;i<fESDOld->GetNumberOfV0s();i++){
1522 AddV0(fESDOld->GetV0(i));
1525 for(int i = 0;i<fESDOld->GetNumberOfCascades();i++){
1526 AddCascade(fESDOld->GetCascade(i));
1529 for(int i = 0;i<fESDOld->GetNumberOfKinks();i++){
1530 AddKink(fESDOld->GetKink(i));
1534 for(int i = 0;i<fESDOld->GetNumberOfCaloClusters();i++){
1535 AddCaloCluster(fESDOld->GetCaloCluster(i));
1541 TObject* AliESDEvent::GetHLTTriggerDecision() const
1543 // get the HLT trigger decission object
1545 // cast away const'nes because the FindListObject method
1547 AliESDEvent* pNonConst=const_cast<AliESDEvent*>(this);
1548 return pNonConst->FindListObject("HLTGlobalTrigger");
1551 TString AliESDEvent::GetHLTTriggerDescription() const
1553 // get the HLT trigger decission description
1554 TString description;
1555 TObject* pDecision=GetHLTTriggerDecision();
1557 description=pDecision->GetTitle();
1563 Bool_t AliESDEvent::IsHLTTriggerFired(const char* name) const
1565 // get the HLT trigger decission description
1566 TObject* pDecision=GetHLTTriggerDecision();
1567 if (!pDecision) return kFALSE;
1569 Option_t* option=pDecision->GetOption();
1570 if (option==NULL || *option!='1') return kFALSE;
1573 TString description=GetHLTTriggerDescription();
1574 Int_t index=description.Index(name);
1575 if (index<0) return kFALSE;
1576 index+=strlen(name);
1577 if (index>=description.Length()) return kFALSE;
1578 if (description[index]!=0 && description[index]!=' ') return kFALSE;