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 <TInterpreter.h>
41 #include "AliESDEvent.h"
42 #include "AliESDfriend.h"
43 #include "AliESDVZERO.h"
44 #include "AliESDFMD.h"
46 #include "AliESDMuonTrack.h"
47 #include "AliESDPmdTrack.h"
48 #include "AliESDTrdTrack.h"
49 #include "AliESDVertex.h"
50 #include "AliESDcascade.h"
51 #include "AliESDPmdTrack.h"
52 #include "AliESDTrdTrigger.h"
53 #include "AliESDTrdTrack.h"
54 #include "AliESDTrdTracklet.h"
55 #include "AliESDVertex.h"
56 #include "AliVertexerTracks.h"
57 #include "AliESDcascade.h"
58 #include "AliESDkink.h"
59 #include "AliESDtrack.h"
60 #include "AliESDHLTtrack.h"
61 #include "AliESDCaloCluster.h"
62 #include "AliESDCaloCells.h"
64 #include "AliESDFMD.h"
65 #include "AliESDVZERO.h"
66 #include "AliMultiplicity.h"
67 #include "AliRawDataErrorLog.h"
69 #include "AliESDACORDE.h"
70 #include "AliESDHLTDecision.h"
71 #include "AliCentrality.h"
72 #include "AliESDCosmicTrack.h"
74 #include "AliESDMFT.h"
76 #include "AliEventplane.h"
83 // here we define the names, some classes are no TNamed, therefore the classnames
85 const char* AliESDEvent::fgkESDListName[kESDListN] = {"AliESDRun",
111 "AliRawDataErrorLogs",
120 //______________________________________________________________________________
121 AliESDEvent::AliESDEvent():
123 fESDObjects(new TList()),
138 fSPDPileupVertices(0),
139 fTrkPileupVertices(0),
149 fEMCALCells(0), fPHOSCells(0),
159 fDetectorStatus(0xFFFFFFFF)
165 //______________________________________________________________________________
166 AliESDEvent::AliESDEvent(const AliESDEvent& esd):
168 fESDObjects(new TList()),
169 fESDRun(new AliESDRun(*esd.fESDRun)),
170 fHeader(new AliESDHeader(*esd.fHeader)),
171 fESDZDC(new AliESDZDC(*esd.fESDZDC)),
172 fESDFMD(new AliESDFMD(*esd.fESDFMD)),
173 fESDVZERO(new AliESDVZERO(*esd.fESDVZERO)),
174 fESDTZERO(new AliESDTZERO(*esd.fESDTZERO)),
175 fTPCVertex(new AliESDVertex(*esd.fTPCVertex)),
176 fSPDVertex(new AliESDVertex(*esd.fSPDVertex)),
177 fPrimaryVertex(new AliESDVertex(*esd.fPrimaryVertex)),
178 fSPDMult(new AliMultiplicity(*esd.fSPDMult)),
179 fPHOSTrigger(new AliESDCaloTrigger(*esd.fPHOSTrigger)),
180 fEMCALTrigger(new AliESDCaloTrigger(*esd.fEMCALTrigger)),
181 fESDACORDE(new AliESDACORDE(*esd.fESDACORDE)),
182 fTrdTrigger(new AliESDTrdTrigger(*esd.fTrdTrigger)),
183 fSPDPileupVertices(new TClonesArray(*esd.fSPDPileupVertices)),
184 fTrkPileupVertices(new TClonesArray(*esd.fTrkPileupVertices)),
185 fTracks(new TClonesArray(*esd.fTracks)),
186 fMuonTracks(new TClonesArray(*esd.fMuonTracks)),
187 fPmdTracks(new TClonesArray(*esd.fPmdTracks)),
188 fTrdTracks(new TClonesArray(*esd.fTrdTracks)),
189 fTrdTracklets(new TClonesArray(*esd.fTrdTracklets)),
190 fV0s(new TClonesArray(*esd.fV0s)),
191 fCascades(new TClonesArray(*esd.fCascades)),
192 fKinks(new TClonesArray(*esd.fKinks)),
193 fCaloClusters(new TClonesArray(*esd.fCaloClusters)),
194 fEMCALCells(new AliESDCaloCells(*esd.fEMCALCells)),
195 fPHOSCells(new AliESDCaloCells(*esd.fPHOSCells)),
196 fCosmicTracks(new TClonesArray(*esd.fCosmicTracks)),
197 fErrorLogs(new TClonesArray(*esd.fErrorLogs)),
198 fESDOld(esd.fESDOld ? new AliESD(*esd.fESDOld) : 0),
199 fESDFriendOld(esd.fESDFriendOld ? new AliESDfriend(*esd.fESDFriendOld) : 0),
200 fConnected(esd.fConnected),
201 fUseOwnList(esd.fUseOwnList),
202 fTOFHeader(new AliTOFHeader(*esd.fTOFHeader)),
203 fCentrality(new AliCentrality(*esd.fCentrality)),
204 fEventplane(new AliEventplane(*esd.fEventplane)),
205 fDetectorStatus(esd.fDetectorStatus)
207 // , fESDMFT(new AliESDMFT(*esd.fESDMFT))
212 printf("copying ESD event...\n"); // AU
213 // CKB init in the constructor list and only add here ...
218 AddObject(fESDVZERO);
219 AddObject(fESDTZERO);
220 AddObject(fTPCVertex);
221 AddObject(fSPDVertex);
222 AddObject(fPrimaryVertex);
224 AddObject(fPHOSTrigger);
225 AddObject(fEMCALTrigger);
226 AddObject(fTrdTrigger);
227 AddObject(fSPDPileupVertices);
228 AddObject(fTrkPileupVertices);
230 AddObject(fMuonTracks);
231 AddObject(fPmdTracks);
232 AddObject(fTrdTracks);
233 AddObject(fTrdTracklets);
235 AddObject(fCascades);
237 AddObject(fCaloClusters);
238 AddObject(fEMCALCells);
239 AddObject(fPHOSCells);
240 AddObject(fCosmicTracks);
241 AddObject(fErrorLogs);
242 AddObject(fESDACORDE);
243 AddObject(fTOFHeader);
245 // AddObject(fESDMFT);
251 //______________________________________________________________________________
252 AliESDEvent & AliESDEvent::operator=(const AliESDEvent& source) {
254 // Assignment operator
256 if(&source == this) return *this;
257 AliVEvent::operator=(source);
259 // This assumes that the list is already created
260 // and that the virtual void Copy(Tobject&) function
261 // is correctly implemented in the derived class
262 // otherwise only TObject::Copy() will be used
266 if((fESDObjects->GetSize()==0)&&(source.fESDObjects->GetSize()>=kESDListN)){
267 // We cover the case that we do not yet have the
268 // standard content but the source has it
272 TIter next(source.GetList());
275 while ((its = next())) {
276 name.Form("%s", its->GetName());
277 TObject *mine = fESDObjects->FindObject(name.Data());
279 TClass* pClass=TClass::GetClass(its->ClassName());
281 AliWarning(Form("Can not find class description for entry %s (%s)\n",
282 its->ClassName(), name.Data()));
286 mine=(TObject*)pClass->New();
288 // not in this: can be added to list
289 AliWarning(Form("%s:%d Could not find %s for copying \n",
290 (char*)__FILE__,__LINE__,name.Data()));
293 if(mine->InheritsFrom("TNamed")){
294 ((TNamed*)mine)->SetName(name);
296 else if(mine->InheritsFrom("TCollection")){
297 if(mine->InheritsFrom("TClonesArray")) {
298 TClonesArray* tcits = dynamic_cast<TClonesArray*>(its);
300 dynamic_cast<TClonesArray*>(mine)->SetClass(tcits->GetClass());
302 dynamic_cast<TCollection*>(mine)->SetName(name);
304 AliDebug(1, Form("adding object %s of type %s", mine->GetName(), mine->ClassName()));
308 if(!its->InheritsFrom("TCollection")){
312 else if(its->InheritsFrom("TClonesArray")){
313 // Create or expand the tclonesarray pointers
314 // so we can directly copy to the object
315 TClonesArray *its_tca = (TClonesArray*)its;
316 TClonesArray *mine_tca = (TClonesArray*)mine;
318 // this leaves the capacity of the TClonesArray the same
319 // except for a factor of 2 increase when size > capacity
320 // does not release any memory occupied by the tca
321 mine_tca->ExpandCreate(its_tca->GetEntriesFast());
322 for(int i = 0;i < its_tca->GetEntriesFast();++i){
324 TObject *mine_tca_obj = mine_tca->At(i);
325 TObject *its_tca_obj = its_tca->At(i);
326 // no need to delete first
327 // pointers within the class should be handled by Copy()...
328 // Can there be Empty slots?
329 its_tca_obj->Copy(*mine_tca_obj);
333 AliWarning(Form("%s:%d cannot copy TCollection \n",
334 (char*)__FILE__,__LINE__));
338 fCentrality = source.fCentrality;
339 fEventplane = source.fEventplane;
341 fConnected = source.fConnected;
342 fUseOwnList = source.fUseOwnList;
344 fDetectorStatus = source.fDetectorStatus;
350 //______________________________________________________________________________
351 AliESDEvent::~AliESDEvent()
354 // Standard destructor
357 // everthing on the list gets deleted automatically
360 if(fESDObjects&&!fConnected)
365 if (fCentrality) delete fCentrality;
366 if (fEventplane) delete fEventplane;
370 void AliESDEvent::Copy(TObject &obj) const {
372 // interface to TOBject::Copy
373 // Copies the content of this into obj!
374 // bascially obj = *this
376 if(this==&obj)return;
377 AliESDEvent *robj = dynamic_cast<AliESDEvent*>(&obj);
378 if(!robj)return; // not an AliESEvent
383 //______________________________________________________________________________
384 void AliESDEvent::Reset()
388 // Std content + Non std content
390 // Reset the standard contents
393 // reset for the old data without AliESDEvent...
394 if(fESDOld)fESDOld->Reset();
396 fESDFriendOld->~AliESDfriend();
397 new (fESDFriendOld) AliESDfriend();
401 if(fESDObjects->GetSize()>kESDListN){
402 // we have non std content
403 // this also covers esdfriends
404 for(int i = kESDListN;i < fESDObjects->GetSize();++i){
405 TObject *pObject = fESDObjects->At(i);
407 if(pObject->InheritsFrom(TClonesArray::Class())){
408 ((TClonesArray*)pObject)->Delete();
410 else if(!pObject->InheritsFrom(TCollection::Class())){
411 TClass *pClass = TClass::GetClass(pObject->ClassName());
412 if (pClass && pClass->GetListOfMethods()->FindObject("Clear")) {
413 AliDebug(1, Form("Clear for object %s class %s", pObject->GetName(), pObject->ClassName()));
417 AliDebug(1, Form("ResetWithPlacementNew for object %s class %s", pObject->GetName(), pObject->ClassName()));
418 ResetWithPlacementNew(pObject);
422 AliWarning(Form("No reset for %s \n",
423 pObject->ClassName()));
430 Bool_t AliESDEvent::ResetWithPlacementNew(TObject *pObject){
431 Long_t dtoronly = TObject::GetDtorOnly();
432 TClass *pClass = TClass::GetClass(pObject->ClassName());
433 TObject::SetDtorOnly(pObject);
435 // Recreate with placement new
436 pClass->New(pObject);
437 // Restore the state.
438 TObject::SetDtorOnly((void*)dtoronly);
442 void AliESDEvent::ResetStdContent()
444 // Reset the standard contents
445 if(fESDRun) fESDRun->Reset();
446 if(fHeader) fHeader->Reset();
447 if(fCentrality) fCentrality->Reset();
448 if(fEventplane) fEventplane->Reset();
449 if(fESDZDC) fESDZDC->Reset();
454 // reset by callin d'to /c'tor keep the pointer
455 fESDVZERO->~AliESDVZERO();
456 new (fESDVZERO) AliESDVZERO();
459 fESDACORDE->~AliESDACORDE();
460 new (fESDACORDE) AliESDACORDE();
462 if(fESDTZERO) fESDTZERO->Reset();
463 // CKB no clear/reset implemented
465 fTPCVertex->~AliESDVertex();
466 new (fTPCVertex) AliESDVertex();
467 fTPCVertex->SetName(fgkESDListName[kTPCVertex]);
470 fSPDVertex->~AliESDVertex();
471 new (fSPDVertex) AliESDVertex();
472 fSPDVertex->SetName(fgkESDListName[kSPDVertex]);
475 fPrimaryVertex->~AliESDVertex();
476 new (fPrimaryVertex) AliESDVertex();
477 fPrimaryVertex->SetName(fgkESDListName[kPrimaryVertex]);
480 fSPDMult->~AliMultiplicity();
481 new (fSPDMult) AliMultiplicity();
484 fTOFHeader->~AliTOFHeader();
485 new (fTOFHeader) AliTOFHeader();
486 //fTOFHeader->SetName(fgkESDListName[kTOFHeader]);
489 fTrdTrigger->~AliESDTrdTrigger();
490 new (fTrdTrigger) AliESDTrdTrigger();
494 // fESDMFT->~AliESDMFT();
495 // new (fESDMFT) AliESDMFT();
499 if(fPHOSTrigger)fPHOSTrigger->DeAllocate();
500 if(fEMCALTrigger)fEMCALTrigger->DeAllocate();
501 if(fSPDPileupVertices)fSPDPileupVertices->Delete();
502 if(fTrkPileupVertices)fTrkPileupVertices->Delete();
503 if(fTracks)fTracks->Delete();
504 if(fMuonTracks)fMuonTracks->Delete();
505 if(fPmdTracks)fPmdTracks->Delete();
506 if(fTrdTracks)fTrdTracks->Delete();
507 if(fTrdTracklets)fTrdTracklets->Delete();
508 if(fV0s)fV0s->Delete();
509 if(fCascades)fCascades->Delete();
510 if(fKinks)fKinks->Delete();
511 if(fCaloClusters)fCaloClusters->Delete();
512 if(fPHOSCells)fPHOSCells->DeleteContainer();
513 if(fEMCALCells)fEMCALCells->DeleteContainer();
514 if(fCosmicTracks)fCosmicTracks->Delete();
515 if(fErrorLogs) fErrorLogs->Delete();
517 // don't reset fconnected fConnected and the list
522 Int_t AliESDEvent::AddV0(const AliESDv0 *v) {
526 TClonesArray &fv = *fV0s;
527 Int_t idx=fV0s->GetEntriesFast();
528 new(fv[idx]) AliESDv0(*v);
532 //______________________________________________________________________________
533 void AliESDEvent::Print(Option_t *) const
536 // Print header information of the event
538 printf("ESD run information\n");
539 printf("Event # in file %d Bunch crossing # %d Orbit # %d Period # %d Run # %d Trigger %lld Magnetic field %f \n",
540 GetEventNumberInFile(),
541 GetBunchCrossNumber(),
546 GetMagneticField() );
548 printf("Vertex: (%.4f +- %.4f, %.4f +- %.4f, %.4f +- %.4f) cm\n",
549 fPrimaryVertex->GetXv(), fPrimaryVertex->GetXRes(),
550 fPrimaryVertex->GetYv(), fPrimaryVertex->GetYRes(),
551 fPrimaryVertex->GetZv(), fPrimaryVertex->GetZRes());
552 printf("Mean vertex in RUN: X=%.4f Y=%.4f Z=%.4f cm\n",
553 GetDiamondX(),GetDiamondY(),GetDiamondZ());
555 printf("SPD Multiplicity. Number of tracklets %d \n",
556 fSPDMult->GetNumberOfTracklets());
557 printf("Number of pileup primary vertices reconstructed with SPD %d\n",
558 GetNumberOfPileupVerticesSPD());
559 printf("Number of pileup primary vertices reconstructed using the tracks %d\n",
560 GetNumberOfPileupVerticesTracks());
561 printf("Number of tracks: \n");
562 printf(" charged %d\n", GetNumberOfTracks());
563 printf(" muon %d\n", GetNumberOfMuonTracks());
564 printf(" pmd %d\n", GetNumberOfPmdTracks());
565 printf(" trd %d\n", GetNumberOfTrdTracks());
566 printf(" trd trkl %d\n", GetNumberOfTrdTracklets());
567 printf(" v0 %d\n", GetNumberOfV0s());
568 printf(" cascades %d\n", GetNumberOfCascades());
569 printf(" kinks %d\n", GetNumberOfKinks());
570 if(fPHOSCells)printf(" PHOSCells %d\n", fPHOSCells->GetNumberOfCells());
571 else printf(" PHOSCells not in the Event\n");
572 if(fEMCALCells)printf(" EMCALCells %d\n", fEMCALCells->GetNumberOfCells());
573 else printf(" EMCALCells not in the Event\n");
574 printf(" CaloClusters %d\n", GetNumberOfCaloClusters());
575 printf(" FMD %s\n", (fESDFMD ? "yes" : "no"));
576 printf(" VZERO %s\n", (fESDVZERO ? "yes" : "no"));
577 if (fCosmicTracks) printf(" Cosmics %d\n", GetNumberOfCosmicTracks());
579 //printf(" MFT %s\n", (fESDMFT ? "yes" : "no"));
583 TObject* pHLTDecision=GetHLTTriggerDecision();
584 printf("HLT trigger decision: %s\n", pHLTDecision?pHLTDecision->GetOption():"not available");
585 if (pHLTDecision) pHLTDecision->Print("compact");
590 void AliESDEvent::SetESDfriend(const AliESDfriend *ev) const {
592 // Attaches the complementary info to the ESD
596 // to be sure that we set the tracks also
597 // in case of old esds
598 // if(fESDOld)CopyFromOldESD();
600 Int_t ntrk=ev->GetNumberOfTracks();
602 for (Int_t i=0; i<ntrk; i++) {
603 const AliESDfriendTrack *f=ev->GetTrack(i);
604 GetTrack(i)->SetFriendTrack(f);
608 Bool_t AliESDEvent::RemoveKink(Int_t rm) const {
609 // ---------------------------------------------------------
610 // Remove a kink candidate and references to it from ESD,
611 // if this candidate does not come from a reconstructed decay
612 // Not yet implemented...
613 // ---------------------------------------------------------
614 Int_t last=GetNumberOfKinks()-1;
615 if ((rm<0)||(rm>last)) return kFALSE;
620 Bool_t AliESDEvent::RemoveV0(Int_t rm) const {
621 // ---------------------------------------------------------
622 // Remove a V0 candidate and references to it from ESD,
623 // if this candidate does not come from a reconstructed decay
624 // ---------------------------------------------------------
625 Int_t last=GetNumberOfV0s()-1;
626 if ((rm<0)||(rm>last)) return kFALSE;
628 AliESDv0 *v0=GetV0(rm);
629 Int_t idxP=v0->GetPindex(), idxN=v0->GetNindex();
632 Int_t lastIdxP=v0->GetPindex(), lastIdxN=v0->GetNindex();
636 // Check if this V0 comes from a reconstructed decay
637 Int_t ncs=GetNumberOfCascades();
638 for (Int_t n=0; n<ncs; n++) {
639 AliESDcascade *cs=GetCascade(n);
641 Int_t csIdxP=cs->GetPindex();
642 Int_t csIdxN=cs->GetNindex();
645 if (idxN==csIdxN) return kFALSE;
647 if (csIdxP==lastIdxP)
648 if (csIdxN==lastIdxN) used++;
651 //Replace the removed V0 with the last V0
652 TClonesArray &a=*fV0s;
653 delete a.RemoveAt(rm);
655 if (rm==last) return kTRUE;
657 //v0 is pointing to the last V0 candidate...
658 new (a[rm]) AliESDv0(*v0);
659 delete a.RemoveAt(last);
661 if (!used) return kTRUE;
664 // Remap the indices of the daughters of reconstructed decays
665 for (Int_t n=0; n<ncs; n++) {
666 AliESDcascade *cs=GetCascade(n);
669 Int_t csIdxP=cs->GetPindex();
670 Int_t csIdxN=cs->GetNindex();
672 if (csIdxP==lastIdxP)
673 if (csIdxN==lastIdxN) {
674 cs->AliESDv0::SetIndex(1,idxP);
675 cs->AliESDv0::SetIndex(0,idxN);
677 if (!used) return kTRUE;
684 Bool_t AliESDEvent::RemoveTrack(Int_t rm) const {
685 // ---------------------------------------------------------
686 // Remove a track and references to it from ESD,
687 // if this track does not come from a reconstructed decay
688 // ---------------------------------------------------------
689 Int_t last=GetNumberOfTracks()-1;
690 if ((rm<0)||(rm>last)) return kFALSE;
694 // Check if this track comes from the reconstructed primary vertices
695 if (fTPCVertex && fTPCVertex->GetStatus()) {
696 UShort_t *primIdx=fTPCVertex->GetIndices();
697 Int_t n=fTPCVertex->GetNIndices();
699 Int_t idx=Int_t(primIdx[n]);
700 if (rm==idx) return kFALSE;
701 if (idx==last) used++;
704 if (fPrimaryVertex && fPrimaryVertex->GetStatus()) {
705 UShort_t *primIdx=fPrimaryVertex->GetIndices();
706 Int_t n=fPrimaryVertex->GetNIndices();
708 Int_t idx=Int_t(primIdx[n]);
709 if (rm==idx) return kFALSE;
710 if (idx==last) used++;
714 // Check if this track comes from a reconstructed decay
715 Int_t nv0=GetNumberOfV0s();
716 for (Int_t n=0; n<nv0; n++) {
717 AliESDv0 *v0=GetV0(n);
719 Int_t idx=v0->GetNindex();
720 if (rm==idx) return kFALSE;
721 if (idx==last) used++;
724 if (rm==idx) return kFALSE;
725 if (idx==last) used++;
728 Int_t ncs=GetNumberOfCascades();
729 for (Int_t n=0; n<ncs; n++) {
730 AliESDcascade *cs=GetCascade(n);
732 Int_t idx=cs->GetIndex();
733 if (rm==idx) return kFALSE;
734 if (idx==last) used++;
738 if (rm==idx) return kFALSE;
739 if (idx==last) used++;
742 if (rm==idx) return kFALSE;
743 if (idx==last) used++;
746 Int_t nkn=GetNumberOfKinks();
747 for (Int_t n=0; n<nkn; n++) {
748 AliESDkink *kn=GetKink(n);
750 Int_t idx=kn->GetIndex(0);
751 if (rm==idx) return kFALSE;
752 if (idx==last) used++;
755 if (rm==idx) return kFALSE;
756 if (idx==last) used++;
759 // Check if this track is associated with a CaloCluster
760 Int_t ncl=GetNumberOfCaloClusters();
761 for (Int_t n=0; n<ncl; n++) {
762 AliESDCaloCluster *cluster=GetCaloCluster(n);
763 TArrayI *arr=cluster->GetTracksMatched();
764 Int_t s=arr->GetSize();
766 Int_t idx=arr->At(s);
767 if (rm==idx) return kFALSE;
768 if (idx==last) used++;
774 //Replace the removed track with the last track
775 TClonesArray &a=*fTracks;
776 delete a.RemoveAt(rm);
778 if (rm==last) return kTRUE;
780 AliESDtrack *t=GetTrack(last);
782 new (a[rm]) AliESDtrack(*t);
783 delete a.RemoveAt(last);
786 if (!used) return kTRUE;
789 // Remap the indices of the tracks used for the primary vertex reconstruction
790 if (fTPCVertex && fTPCVertex->GetStatus()) {
791 UShort_t *primIdx=fTPCVertex->GetIndices();
792 Int_t n=fTPCVertex->GetNIndices();
794 Int_t idx=Int_t(primIdx[n]);
796 primIdx[n]=Short_t(rm);
798 if (!used) return kTRUE;
802 if (fPrimaryVertex && fPrimaryVertex->GetStatus()) {
803 UShort_t *primIdx=fPrimaryVertex->GetIndices();
804 Int_t n=fPrimaryVertex->GetNIndices();
806 Int_t idx=Int_t(primIdx[n]);
808 primIdx[n]=Short_t(rm);
810 if (!used) return kTRUE;
815 // Remap the indices of the daughters of reconstructed decays
816 for (Int_t n=0; n<nv0; n++) {
817 AliESDv0 *v0=GetV0(n);
818 if (v0->GetIndex(0)==last) {
821 if (!used) return kTRUE;
823 if (v0->GetIndex(1)==last) {
826 if (!used) return kTRUE;
830 for (Int_t n=0; n<ncs; n++) {
831 AliESDcascade *cs=GetCascade(n);
832 if (cs->GetIndex()==last) {
835 if (!used) return kTRUE;
838 if (v0->GetIndex(0)==last) {
841 if (!used) return kTRUE;
843 if (v0->GetIndex(1)==last) {
846 if (!used) return kTRUE;
850 for (Int_t n=0; n<nkn; n++) {
851 AliESDkink *kn=GetKink(n);
852 if (kn->GetIndex(0)==last) {
855 if (!used) return kTRUE;
857 if (kn->GetIndex(1)==last) {
860 if (!used) return kTRUE;
864 // Remap the indices of the tracks accosicated with CaloClusters
865 for (Int_t n=0; n<ncl; n++) {
866 AliESDCaloCluster *cluster=GetCaloCluster(n);
867 TArrayI *arr=cluster->GetTracksMatched();
868 Int_t s=arr->GetSize();
870 Int_t idx=arr->At(s);
874 if (!used) return kTRUE;
883 Bool_t AliESDEvent::Clean(Float_t *cleanPars) {
885 // Remove the data which are not needed for the physics analysis.
887 // 1) Cleaning the V0 candidates
888 // ---------------------------
889 // If the cosine of the V0 pointing angle "csp" and
890 // the DCA between the daughter tracks "dca" does not satisfy
893 // csp > cleanPars[1] + dca/cleanPars[0]*(1.- cleanPars[1])
895 // an attempt to remove this V0 candidate from ESD is made.
897 // The V0 candidate gets removed if it does not belong to any
898 // recosntructed cascade decay
900 // 12.11.2007, optimal values: cleanPars[0]=0.5, cleanPars[1]=0.999
902 // 2) Cleaning the tracks
903 // ----------------------
904 // If track's transverse parameter is larger than cleanPars[2]
906 // track's longitudinal parameter is larger than cleanPars[3]
907 // an attempt to remove this track from ESD is made.
909 // The track gets removed if it does not come
910 // from a reconstructed decay
914 Float_t dcaMax=cleanPars[0];
915 Float_t cspMin=cleanPars[1];
917 Int_t nV0s=GetNumberOfV0s();
918 for (Int_t i=nV0s-1; i>=0; i--) {
919 AliESDv0 *v0=GetV0(i);
921 Float_t dca=v0->GetDcaV0Daughters();
922 Float_t csp=v0->GetV0CosineOfPointingAngle();
923 Float_t cspcut=cspMin + dca/dcaMax*(1.-cspMin);
924 if (csp > cspcut) continue;
925 if (RemoveV0(i)) rc=kTRUE;
929 Float_t dmax=cleanPars[2], zmax=cleanPars[3];
931 const AliESDVertex *vertex=GetPrimaryVertexSPD();
932 Bool_t vtxOK=vertex->GetStatus();
934 Int_t nTracks=GetNumberOfTracks();
935 for (Int_t i=nTracks-1; i>=0; i--) {
936 AliESDtrack *track=GetTrack(i);
937 Float_t xy,z; track->GetImpactParameters(xy,z);
938 if ((TMath::Abs(xy) > dmax) || (vtxOK && (TMath::Abs(z) > zmax))) {
939 if (RemoveTrack(i)) rc=kTRUE;
946 Char_t AliESDEvent::AddPileupVertexSPD(const AliESDVertex *vtx)
948 // Add a pileup primary vertex reconstructed with SPD
949 TClonesArray &ftr = *fSPDPileupVertices;
950 Char_t n=Char_t(ftr.GetEntriesFast());
951 AliESDVertex *vertex = new(ftr[n]) AliESDVertex(*vtx);
956 Char_t AliESDEvent::AddPileupVertexTracks(const AliESDVertex *vtx)
958 // Add a pileup primary vertex reconstructed with SPD
959 TClonesArray &ftr = *fTrkPileupVertices;
960 Char_t n=Char_t(ftr.GetEntriesFast());
961 AliESDVertex *vertex = new(ftr[n]) AliESDVertex(*vtx);
966 Int_t AliESDEvent::AddTrack(const AliESDtrack *t)
969 TClonesArray &ftr = *fTracks;
970 AliESDtrack * track = new(ftr[fTracks->GetEntriesFast()])AliESDtrack(*t);
971 track->SetID(fTracks->GetEntriesFast()-1);
972 return track->GetID();
975 AliESDtrack* AliESDEvent::NewTrack()
978 TClonesArray &ftr = *fTracks;
979 AliESDtrack * track = new(ftr[fTracks->GetEntriesFast()])AliESDtrack();
980 track->SetID(fTracks->GetEntriesFast()-1);
984 void AliESDEvent::AddMuonTrack(const AliESDMuonTrack *t)
986 TClonesArray &fmu = *fMuonTracks;
987 new(fmu[fMuonTracks->GetEntriesFast()]) AliESDMuonTrack(*t);
990 void AliESDEvent::AddPmdTrack(const AliESDPmdTrack *t)
992 TClonesArray &fpmd = *fPmdTracks;
993 new(fpmd[fPmdTracks->GetEntriesFast()]) AliESDPmdTrack(*t);
996 void AliESDEvent::SetTrdTrigger(const AliESDTrdTrigger *t)
1001 void AliESDEvent::AddTrdTrack(const AliESDTrdTrack *t)
1003 TClonesArray &ftrd = *fTrdTracks;
1004 new(ftrd[fTrdTracks->GetEntriesFast()]) AliESDTrdTrack(*t);
1007 void AliESDEvent::AddTrdTracklet(const AliESDTrdTracklet *trkl)
1009 new ((*fTrdTracklets)[fTrdTracklets->GetEntriesFast()]) AliESDTrdTracklet(*trkl);
1012 Int_t AliESDEvent::AddKink(const AliESDkink *c)
1015 TClonesArray &fk = *fKinks;
1016 AliESDkink * kink = new(fk[fKinks->GetEntriesFast()]) AliESDkink(*c);
1017 kink->SetID(fKinks->GetEntriesFast()); // CKB different from the other imps..
1018 return fKinks->GetEntriesFast()-1;
1022 void AliESDEvent::AddCascade(const AliESDcascade *c)
1024 TClonesArray &fc = *fCascades;
1025 new(fc[fCascades->GetEntriesFast()]) AliESDcascade(*c);
1028 void AliESDEvent::AddCosmicTrack(const AliESDCosmicTrack *t)
1030 TClonesArray &ft = *fCosmicTracks;
1031 new(ft[fCosmicTracks->GetEntriesFast()]) AliESDCosmicTrack(*t);
1035 Int_t AliESDEvent::AddCaloCluster(const AliESDCaloCluster *c)
1038 TClonesArray &fc = *fCaloClusters;
1039 AliESDCaloCluster *clus = new(fc[fCaloClusters->GetEntriesFast()]) AliESDCaloCluster(*c);
1040 clus->SetID(fCaloClusters->GetEntriesFast()-1);
1041 return fCaloClusters->GetEntriesFast()-1;
1045 void AliESDEvent::AddRawDataErrorLog(const AliRawDataErrorLog *log) const {
1046 TClonesArray &errlogs = *fErrorLogs;
1047 new(errlogs[errlogs.GetEntriesFast()]) AliRawDataErrorLog(*log);
1050 void AliESDEvent::SetZDCData(AliESDZDC * obj)
1052 // use already allocated space
1057 void AliESDEvent::SetPrimaryVertexTPC(const AliESDVertex *vertex)
1059 // Set the TPC vertex
1060 // use already allocated space
1062 *fTPCVertex = *vertex;
1063 fTPCVertex->SetName(fgkESDListName[kTPCVertex]);
1067 void AliESDEvent::SetPrimaryVertexSPD(const AliESDVertex *vertex)
1069 // Set the SPD vertex
1070 // use already allocated space
1072 *fSPDVertex = *vertex;
1073 fSPDVertex->SetName(fgkESDListName[kSPDVertex]);
1077 void AliESDEvent::SetPrimaryVertexTracks(const AliESDVertex *vertex)
1079 // Set the primary vertex reconstructed using he ESD tracks.
1080 // use already allocated space
1082 *fPrimaryVertex = *vertex;
1083 fPrimaryVertex->SetName(fgkESDListName[kPrimaryVertex]);
1087 const AliESDVertex * AliESDEvent::GetPrimaryVertex() const
1090 // Get the "best" available reconstructed primary vertex.
1093 if (fPrimaryVertex->GetStatus()) return fPrimaryVertex;
1096 if (fSPDVertex->GetStatus()) return fSPDVertex;
1098 if(fTPCVertex) return fTPCVertex;
1100 AliWarning("No primary vertex available. Returning the \"default\"...");
1104 AliESDVertex * AliESDEvent::PrimaryVertexTracksUnconstrained() const
1107 // Removes diamond constraint from fPrimaryVertex (reconstructed with tracks)
1108 // Returns a AliESDVertex which has to be deleted by the user
1110 if(!fPrimaryVertex) {
1111 AliWarning("No primary vertex from tracks available.");
1114 if(!fPrimaryVertex->GetStatus()) {
1115 AliWarning("No primary vertex from tracks available.");
1119 AliVertexerTracks vertexer(GetMagneticField());
1120 Float_t diamondxyz[3]={(Float_t)GetDiamondX(),(Float_t)GetDiamondY(),0.};
1121 Float_t diamondcovxy[3]; GetDiamondCovXY(diamondcovxy);
1122 Float_t diamondcov[6]={diamondcovxy[0],diamondcovxy[1],diamondcovxy[2],0.,0.,7.};
1123 AliESDVertex *vertex =
1124 (AliESDVertex*)vertexer.RemoveConstraintFromVertex(fPrimaryVertex,diamondxyz,diamondcov);
1129 void AliESDEvent::SetMultiplicity(const AliMultiplicity *mul)
1131 // Set the SPD Multiplicity
1138 void AliESDEvent::SetFMDData(AliESDFMD * obj)
1140 // use already allocated space
1146 void AliESDEvent::SetVZEROData(AliESDVZERO * obj)
1148 // use already allocated space
1153 void AliESDEvent::SetTZEROData(AliESDTZERO * obj)
1155 // use already allocated space
1161 //void AliESDEvent::SetMFTData(AliESDMFT * obj)
1168 void AliESDEvent::SetACORDEData(AliESDACORDE * obj)
1175 void AliESDEvent::GetESDfriend(AliESDfriend *ev) const
1178 // Extracts the complementary info from the ESD
1182 Int_t ntrk=GetNumberOfTracks();
1184 for (Int_t i=0; i<ntrk; i++) {
1185 AliESDtrack *t=GetTrack(i);
1186 const AliESDfriendTrack *f=t->GetFriendTrack();
1189 t->ReleaseESDfriendTrack();// Not to have two copies of "friendTrack"
1193 AliESDfriend *fr = (AliESDfriend*)(const_cast<AliESDEvent*>(this)->FindListObject("AliESDfriend"));
1194 if (fr) ev->SetVZEROfriend(fr->GetVZEROfriend());
1197 void AliESDEvent::AddObject(TObject* obj)
1199 // Add an object to the list of object.
1200 // Please be aware that in order to increase performance you should
1201 // refrain from using TObjArrays (if possible). Use TClonesArrays, instead.
1202 fESDObjects->SetOwner(kTRUE);
1203 fESDObjects->AddLast(obj);
1207 void AliESDEvent::GetStdContent()
1209 // set pointers for standard content
1210 // get by name much safer and not a big overhead since not called very often
1212 fESDRun = (AliESDRun*)fESDObjects->FindObject(fgkESDListName[kESDRun]);
1213 fHeader = (AliESDHeader*)fESDObjects->FindObject(fgkESDListName[kHeader]);
1214 fESDZDC = (AliESDZDC*)fESDObjects->FindObject(fgkESDListName[kESDZDC]);
1215 fESDFMD = (AliESDFMD*)fESDObjects->FindObject(fgkESDListName[kESDFMD]);
1216 fESDVZERO = (AliESDVZERO*)fESDObjects->FindObject(fgkESDListName[kESDVZERO]);
1217 fESDTZERO = (AliESDTZERO*)fESDObjects->FindObject(fgkESDListName[kESDTZERO]);
1218 fTPCVertex = (AliESDVertex*)fESDObjects->FindObject(fgkESDListName[kTPCVertex]);
1219 fSPDVertex = (AliESDVertex*)fESDObjects->FindObject(fgkESDListName[kSPDVertex]);
1220 fPrimaryVertex = (AliESDVertex*)fESDObjects->FindObject(fgkESDListName[kPrimaryVertex]);
1221 fSPDMult = (AliMultiplicity*)fESDObjects->FindObject(fgkESDListName[kSPDMult]);
1222 fPHOSTrigger = (AliESDCaloTrigger*)fESDObjects->FindObject(fgkESDListName[kPHOSTrigger]);
1223 fEMCALTrigger = (AliESDCaloTrigger*)fESDObjects->FindObject(fgkESDListName[kEMCALTrigger]);
1224 fSPDPileupVertices = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kSPDPileupVertices]);
1225 fTrkPileupVertices = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kTrkPileupVertices]);
1226 fTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kTracks]);
1227 fMuonTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kMuonTracks]);
1228 fPmdTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kPmdTracks]);
1229 fTrdTrigger = (AliESDTrdTrigger*)fESDObjects->FindObject(fgkESDListName[kTrdTrigger]);
1230 fTrdTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kTrdTracks]);
1231 fTrdTracklets = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kTrdTracklets]);
1232 fV0s = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kV0s]);
1233 fCascades = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kCascades]);
1234 fKinks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kKinks]);
1235 fCaloClusters = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kCaloClusters]);
1236 fEMCALCells = (AliESDCaloCells*)fESDObjects->FindObject(fgkESDListName[kEMCALCells]);
1237 fPHOSCells = (AliESDCaloCells*)fESDObjects->FindObject(fgkESDListName[kPHOSCells]);
1238 fErrorLogs = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kErrorLogs]);
1239 fESDACORDE = (AliESDACORDE*)fESDObjects->FindObject(fgkESDListName[kESDACORDE]);
1240 fTOFHeader = (AliTOFHeader*)fESDObjects->FindObject(fgkESDListName[kTOFHeader]);
1241 fCosmicTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kCosmicTracks]);
1243 // fESDMFT = (AliESDMFT*)fESDObjects->FindObject(fgkESDListName[kESDMFT]);
1247 void AliESDEvent::SetStdNames(){
1248 // Set the names of the standard contents
1250 if(fESDObjects->GetEntries()>=kESDListN){
1251 for(int i = 0;i < fESDObjects->GetEntries() && i<kESDListN;i++){
1252 TObject *fObj = fESDObjects->At(i);
1253 if(fObj->InheritsFrom("TNamed")){
1254 ((TNamed*)fObj)->SetName(fgkESDListName[i]);
1256 else if(fObj->InheritsFrom("TClonesArray")){
1257 ((TClonesArray*)fObj)->SetName(fgkESDListName[i]);
1262 AliWarning("Std Entries missing");
1267 void AliESDEvent::CreateStdContent(Bool_t bUseThisList){
1268 fUseOwnList = bUseThisList;
1272 void AliESDEvent::CreateStdContent()
1274 // create the standard AOD content and set pointers
1276 // create standard objects and add them to the TList of objects
1277 AddObject(new AliESDRun());
1278 AddObject(new AliESDHeader());
1279 AddObject(new AliESDZDC());
1280 AddObject(new AliESDFMD());
1281 AddObject(new AliESDVZERO());
1282 AddObject(new AliESDTZERO());
1283 AddObject(new AliESDVertex());
1284 AddObject(new AliESDVertex());
1285 AddObject(new AliESDVertex());
1286 AddObject(new AliMultiplicity());
1287 AddObject(new AliESDCaloTrigger());
1288 AddObject(new AliESDCaloTrigger());
1289 AddObject(new TClonesArray("AliESDVertex",0));
1290 AddObject(new TClonesArray("AliESDVertex",0));
1291 AddObject(new TClonesArray("AliESDtrack",0));
1292 AddObject(new TClonesArray("AliESDMuonTrack",0));
1293 AddObject(new TClonesArray("AliESDPmdTrack",0));
1294 AddObject(new AliESDTrdTrigger());
1295 AddObject(new TClonesArray("AliESDTrdTrack",0));
1296 AddObject(new TClonesArray("AliESDTrdTracklet",0));
1297 AddObject(new TClonesArray("AliESDv0",0));
1298 AddObject(new TClonesArray("AliESDcascade",0));
1299 AddObject(new TClonesArray("AliESDkink",0));
1300 AddObject(new TClonesArray("AliESDCaloCluster",0));
1301 AddObject(new AliESDCaloCells());
1302 AddObject(new AliESDCaloCells());
1303 AddObject(new TClonesArray("AliRawDataErrorLog",0));
1304 AddObject(new AliESDACORDE());
1305 AddObject(new AliTOFHeader());
1306 AddObject(new TClonesArray("AliESDCosmicTrack",0));
1308 //AddObject(new AliESDMFT());
1311 // check the order of the indices against enum...
1315 // read back pointers
1319 TObject* AliESDEvent::FindListObject(const char *name) const {
1321 // Find object with name "name" in the list of branches
1324 return fESDObjects->FindObject(name);
1329 Int_t AliESDEvent::GetPHOSClusters(TRefArray *clusters) const
1331 // fills the provided TRefArray with all found phos clusters
1335 AliESDCaloCluster *cl = 0;
1336 for (Int_t i = 0; i < GetNumberOfCaloClusters(); i++) {
1338 if ( (cl = GetCaloCluster(i)) ) {
1341 AliDebug(1,Form("IsPHOS cluster %d Size: %d \n",i,clusters->GetEntriesFast()));
1345 return clusters->GetEntriesFast();
1348 Int_t AliESDEvent::GetEMCALClusters(TRefArray *clusters) const
1350 // fills the provided TRefArray with all found emcal clusters
1354 AliESDCaloCluster *cl = 0;
1355 for (Int_t i = 0; i < GetNumberOfCaloClusters(); i++) {
1357 if ( (cl = GetCaloCluster(i)) ) {
1360 AliDebug(1,Form("IsEMCAL cluster %d Size: %d \n",i,clusters->GetEntriesFast()));
1364 return clusters->GetEntriesFast();
1367 void AliESDEvent::WriteToTree(TTree* tree) const {
1368 // Book the branches as in TTree::Branch(TCollection*)
1369 // but add a "." at the end of top level branches which are
1370 // not a TClonesArray
1374 TIter next(fESDObjects);
1375 const Int_t kSplitlevel = 99; // default value in TTree::Branch()
1376 const Int_t kBufsize = 32000; // default value in TTree::Branch()
1379 while ((obj = next())) {
1380 branchname.Form("%s", obj->GetName());
1381 if(branchname.CompareTo("AliESDfriend")==0)branchname = "ESDfriend.";
1382 if ((kSplitlevel > 1) && !obj->InheritsFrom(TClonesArray::Class())) {
1383 if(!branchname.EndsWith("."))branchname += ".";
1385 if (!tree->FindBranch(branchname)) {
1386 // For the custom streamer to be called splitlevel
1387 // has to be negative, only needed for HLT
1388 Int_t splitLevel = (TString(obj->ClassName()) == "AliHLTGlobalTriggerDecision") ? -1 : kSplitlevel - 1;
1389 tree->Bronch(branchname, obj->ClassName(), fESDObjects->GetObjectRef(obj),kBufsize, splitLevel);
1395 void AliESDEvent::ReadFromTree(TTree *tree, Option_t* opt){
1397 // Connect the ESDEvent to a tree
1400 AliWarning("AliESDEvent::ReadFromTree() Zero Pointer to Tree \n");
1404 if(!tree->GetTree())tree->LoadTree(0);
1406 // if we find the "ESD" branch on the tree we do have the old structure
1407 if(tree->GetBranch("ESD")) {
1408 char ** address = (char **)(tree->GetBranch("ESD")->GetAddress());
1409 // do we have the friend branch
1410 TBranch * esdFB = tree->GetBranch("ESDfriend.");
1411 char ** addressF = 0;
1412 if(esdFB)addressF = (char **)(esdFB->GetAddress());
1414 AliInfo("AliESDEvent::ReadFromTree() Reading old Tree");
1415 tree->SetBranchAddress("ESD", &fESDOld);
1417 tree->SetBranchAddress("ESDfriend.",&fESDFriendOld);
1420 AliInfo("AliESDEvent::ReadFromTree() Reading old Tree");
1421 AliInfo("Branch already connected. Using existing branch address.");
1422 fESDOld = (AliESD*) (*address);
1423 // addressF can still be 0, since branch needs to switched on
1424 if(addressF)fESDFriendOld = (AliESDfriend*) (*addressF);
1427 // have already connected the old ESD structure... ?
1428 // reuse also the pointer of the AlliESDEvent
1429 // otherwise create new ones
1430 TList* connectedList = (TList*) (tree->GetUserInfo()->FindObject("ESDObjectsConnectedToTree"));
1433 // If connected use the connected list of objects
1434 if(fESDObjects!= connectedList){
1435 // protect when called twice
1436 fESDObjects->Delete();
1437 fESDObjects = connectedList;
1442 // The pointer to the friend changes when called twice via InitIO
1443 // since AliESDEvent is deleted
1444 TObject* oldf = FindListObject("AliESDfriend");
1447 newf = (TObject*)*addressF;
1449 if(newf!=0&&oldf!=newf){
1450 // remove the old reference
1451 // Should we also delete it? Or is this handled in TTree I/O
1452 // since it is created by the first SetBranchAddress
1453 fESDObjects->Remove(oldf);
1455 fESDObjects->Add(newf);
1462 CreateStdContent(); // create for copy
1463 // if we have the esdfriend add it, so we always can access it via the userinfo
1464 if(fESDFriendOld)AddObject(fESDFriendOld);
1465 // we are not owner of the list objects
1466 // must not delete it
1467 fESDObjects->SetOwner(kTRUE);
1468 fESDObjects->SetName("ESDObjectsConnectedToTree");
1469 tree->GetUserInfo()->Add(fESDObjects);
1477 // Try to find AliESDEvent
1478 AliESDEvent *esdEvent = 0;
1479 esdEvent = (AliESDEvent*)tree->GetTree()->GetUserInfo()->FindObject("AliESDEvent");
1481 // Check if already connected to tree
1483 TList* connectedList = (TList*) (tree->GetUserInfo()->FindObject("ESDObjectsConnectedToTree"));
1486 if (connectedList && (strcmp(opt, "reconnect"))) {
1487 // If connected use the connected list if objects
1488 fESDObjects->Delete();
1489 fESDObjects = connectedList;
1496 // prevent a memory leak when reading back the TList
1497 // if (!(strcmp(opt, "reconnect"))) fESDObjects->Delete();
1500 // create a new TList from the UserInfo TList...
1501 // copy constructor does not work...
1502 fESDObjects = (TList*)(esdEvent->GetList()->Clone());
1503 fESDObjects->SetOwner(kTRUE);
1505 else if ( fESDObjects->GetEntries()==0){
1506 // at least create the std content if we want to read to our list
1511 // we only need new things in the list if we do no already have it..
1512 // TODO just add new entries
1514 if(fESDObjects->GetEntries()<kESDListN){
1515 AliWarning(Form("AliESDEvent::ReadFromTree() TList contains less than the standard contents %d < %d \n",
1516 fESDObjects->GetEntries(),kESDListN));
1518 // set the branch addresses
1519 TIter next(fESDObjects);
1521 while((el=(TNamed*)next())){
1522 TString bname(el->GetName());
1523 if(bname.CompareTo("AliESDfriend")==0)
1525 // AliESDfriend does not have a name ...
1526 TBranch *br = tree->GetBranch("ESDfriend.");
1527 if (br) tree->SetBranchAddress("ESDfriend.",fESDObjects->GetObjectRef(el));
1530 // check if branch exists under this Name
1531 TBranch *br = tree->GetBranch(bname.Data());
1533 tree->SetBranchAddress(bname.Data(),fESDObjects->GetObjectRef(el));
1536 br = tree->GetBranch(Form("%s.",bname.Data()));
1538 tree->SetBranchAddress(Form("%s.",bname.Data()),fESDObjects->GetObjectRef(el));
1541 AliWarning(Form("AliESDEvent::ReadFromTree() No Branch found with Name %s or %s.",bname.Data(),bname.Data()));
1548 // when reading back we are not owner of the list
1549 // must not delete it
1550 fESDObjects->SetOwner(kTRUE);
1551 fESDObjects->SetName("ESDObjectsConnectedToTree");
1552 // we are not owner of the list objects
1553 // must not delete it
1554 tree->GetUserInfo()->Add(fESDObjects);
1555 tree->GetUserInfo()->SetOwner(kFALSE);
1559 // we can't get the list from the user data, create standard content
1560 // and set it by hand (no ESDfriend at the moment
1562 TIter next(fESDObjects);
1564 while((el=(TNamed*)next())){
1565 TString bname(el->GetName());
1566 TBranch *br = tree->GetBranch(bname.Data());
1568 tree->SetBranchAddress(bname.Data(),fESDObjects->GetObjectRef(el));
1571 br = tree->GetBranch(Form("%s.",bname.Data()));
1573 tree->SetBranchAddress(Form("%s.",bname.Data()),fESDObjects->GetObjectRef(el));
1578 // when reading back we are not owner of the list
1579 // must not delete it
1580 fESDObjects->SetOwner(kTRUE);
1585 void AliESDEvent::CopyFromOldESD()
1587 // Method which copies over everthing from the old esd structure to the
1592 SetRunNumber(fESDOld->GetRunNumber());
1593 SetPeriodNumber(fESDOld->GetPeriodNumber());
1594 SetMagneticField(fESDOld->GetMagneticField());
1596 // leave out diamond ...
1597 // SetDiamond(const AliESDVertex *vertex) { fESDRun->SetDiamond(vertex);}
1600 SetTriggerMask(fESDOld->GetTriggerMask());
1601 SetOrbitNumber(fESDOld->GetOrbitNumber());
1602 SetTimeStamp(fESDOld->GetTimeStamp());
1603 SetEventType(fESDOld->GetEventType());
1604 SetEventNumberInFile(fESDOld->GetEventNumberInFile());
1605 SetBunchCrossNumber(fESDOld->GetBunchCrossNumber());
1606 SetTriggerCluster(fESDOld->GetTriggerCluster());
1610 SetZDC(fESDOld->GetZDCN1Energy(),
1611 fESDOld->GetZDCP1Energy(),
1612 fESDOld->GetZDCEMEnergy(),
1614 fESDOld->GetZDCN2Energy(),
1615 fESDOld->GetZDCP2Energy(),
1616 fESDOld->GetZDCParticipants(),
1626 if(fESDOld->GetFMDData())SetFMDData(fESDOld->GetFMDData());
1630 SetT0zVertex(fESDOld->GetT0zVertex());
1631 SetT0(fESDOld->GetT0());
1635 if (fESDOld->GetVZEROData()) SetVZEROData(fESDOld->GetVZEROData());
1637 if(fESDOld->GetVertex())SetPrimaryVertexSPD(fESDOld->GetVertex());
1639 if(fESDOld->GetPrimaryVertex())SetPrimaryVertexTracks(fESDOld->GetPrimaryVertex());
1641 if(fESDOld->GetMultiplicity())SetMultiplicity(fESDOld->GetMultiplicity());
1643 for(int i = 0;i<fESDOld->GetNumberOfTracks();i++){
1644 AddTrack(fESDOld->GetTrack(i));
1647 for(int i = 0;i<fESDOld->GetNumberOfMuonTracks();i++){
1648 AddMuonTrack(fESDOld->GetMuonTrack(i));
1651 for(int i = 0;i<fESDOld->GetNumberOfPmdTracks();i++){
1652 AddPmdTrack(fESDOld->GetPmdTrack(i));
1655 for(int i = 0;i<fESDOld->GetNumberOfTrdTracks();i++){
1656 AddTrdTrack(fESDOld->GetTrdTrack(i));
1659 for(int i = 0;i<fESDOld->GetNumberOfV0s();i++){
1660 AddV0(fESDOld->GetV0(i));
1663 for(int i = 0;i<fESDOld->GetNumberOfCascades();i++){
1664 AddCascade(fESDOld->GetCascade(i));
1667 for(int i = 0;i<fESDOld->GetNumberOfKinks();i++){
1668 AddKink(fESDOld->GetKink(i));
1672 for(int i = 0;i<fESDOld->GetNumberOfCaloClusters();i++){
1673 AddCaloCluster(fESDOld->GetCaloCluster(i));
1678 // if (fESDOld->GetMFTData()) SetMFTData(fESDOld->GetMFTData());
1684 Bool_t AliESDEvent::IsEventSelected(const char *trigExpr) const
1686 // Check if the event satisfies the trigger
1687 // selection expression trigExpr.
1688 // trigExpr can be any logical expression
1689 // of the trigger classes defined in AliESDRun
1690 // In case of wrong syntax return kTRUE.
1692 TString expr(trigExpr);
1693 if (expr.IsNull()) return kTRUE;
1695 ULong64_t mask = GetTriggerMask();
1696 for(Int_t itrig = 0; itrig < AliESDRun::kNTriggerClasses; itrig++) {
1697 if (mask & (1ull << itrig)) {
1698 expr.ReplaceAll(GetESDRun()->GetTriggerClass(itrig),"1");
1701 expr.ReplaceAll(GetESDRun()->GetTriggerClass(itrig),"0");
1706 if ((gROOT->ProcessLineFast(expr.Data(),&error) == 0) &&
1707 (error == TInterpreter::kNoError)) {
1715 TObject* AliESDEvent::GetHLTTriggerDecision() const
1717 // get the HLT trigger decission object
1719 // cast away const'nes because the FindListObject method
1721 AliESDEvent* pNonConst=const_cast<AliESDEvent*>(this);
1722 return pNonConst->FindListObject("HLTGlobalTrigger");
1725 TString AliESDEvent::GetHLTTriggerDescription() const
1727 // get the HLT trigger decission description
1728 TString description;
1729 TObject* pDecision=GetHLTTriggerDecision();
1731 description=pDecision->GetTitle();
1737 Bool_t AliESDEvent::IsHLTTriggerFired(const char* name) const
1739 // get the HLT trigger decission description
1740 TObject* pDecision=GetHLTTriggerDecision();
1741 if (!pDecision) return kFALSE;
1743 Option_t* option=pDecision->GetOption();
1744 if (option==NULL || *option!='1') return kFALSE;
1747 TString description=GetHLTTriggerDescription();
1748 Int_t index=description.Index(name);
1749 if (index<0) return kFALSE;
1750 index+=strlen(name);
1751 if (index>=description.Length()) return kFALSE;
1752 if (description[index]!=0 && description[index]!=' ') return kFALSE;
1757 //______________________________________________________________________________
1758 Bool_t AliESDEvent::IsPileupFromSPD(Int_t minContributors,
1760 Double_t nSigmaZdist,
1761 Double_t nSigmaDiamXY,
1762 Double_t nSigmaDiamZ) const{
1764 // This function checks if there was a pile up
1765 // reconstructed with SPD
1767 Int_t nc1=fSPDVertex->GetNContributors();
1768 if(nc1<1) return kFALSE;
1769 Int_t nPileVert=GetNumberOfPileupVerticesSPD();
1770 if(nPileVert==0) return kFALSE;
1772 for(Int_t i=0; i<nPileVert;i++){
1773 const AliESDVertex* pv=GetPileupVertexSPD(i);
1774 Int_t nc2=pv->GetNContributors();
1775 if(nc2>=minContributors){
1776 Double_t z1=fSPDVertex->GetZ();
1777 Double_t z2=pv->GetZ();
1778 Double_t distZ=TMath::Abs(z2-z1);
1779 Double_t distZdiam=TMath::Abs(z2-GetDiamondZ());
1780 Double_t cutZdiam=nSigmaDiamZ*TMath::Sqrt(GetSigma2DiamondZ());
1781 if(GetSigma2DiamondZ()<0.0001)cutZdiam=99999.; //protection for missing z diamond information
1782 if(distZ>minZdist && distZdiam<cutZdiam){
1783 Double_t x2=pv->GetX();
1784 Double_t y2=pv->GetY();
1785 Double_t distXdiam=TMath::Abs(x2-GetDiamondX());
1786 Double_t distYdiam=TMath::Abs(y2-GetDiamondY());
1787 Double_t cov1[6],cov2[6];
1788 fSPDVertex->GetCovarianceMatrix(cov1);
1789 pv->GetCovarianceMatrix(cov2);
1790 Double_t errxDist=TMath::Sqrt(cov2[0]+GetSigma2DiamondX());
1791 Double_t erryDist=TMath::Sqrt(cov2[2]+GetSigma2DiamondY());
1792 Double_t errzDist=TMath::Sqrt(cov1[5]+cov2[5]);
1793 Double_t cutXdiam=nSigmaDiamXY*errxDist;
1794 if(GetSigma2DiamondX()<0.0001)cutXdiam=99999.; //protection for missing diamond information
1795 Double_t cutYdiam=nSigmaDiamXY*erryDist;
1796 if(GetSigma2DiamondY()<0.0001)cutYdiam=99999.; //protection for missing diamond information
1797 if( (distXdiam<cutXdiam) && (distYdiam<cutYdiam) && (distZ>nSigmaZdist*errzDist) ){
1806 //______________________________________________________________________________
1807 void AliESDEvent::EstimateMultiplicity(Int_t &tracklets, Int_t &trITSTPC, Int_t &trITSSApure, Double_t eta, Bool_t useDCAFlag,Bool_t useV0Flag) const
1810 // calculates 3 estimators for the multiplicity in the -eta:eta range
1811 // tracklets : using SPD tracklets only
1812 // trITSTPC : using TPC/ITS + complementary ITS SA tracks + tracklets from clusters not used by tracks
1813 // trITSSApure : using ITS standalone tracks + tracklets from clusters not used by tracks
1814 // if useDCAFlag is true: account for the ESDtrack flag marking the tracks with large DCA
1815 // if useV0Flag is true: account for the ESDtrack flag marking conversion and K0's V0s
1816 tracklets = trITSSApure = trITSTPC = 0;
1817 int ntr = fSPDMult ? fSPDMult->GetNumberOfTracklets() : 0;
1820 for (int itr=ntr;itr--;) {
1821 if (TMath::Abs(fSPDMult->GetEta(itr))>eta) continue;
1823 if (fSPDMult->FreeClustersTracklet(itr,0)) trITSTPC++; // not used in ITS/TPC or ITS_SA track
1824 if (fSPDMult->FreeClustersTracklet(itr,1)) trITSSApure++; // not used in ITS_SA_Pure track
1827 // count real tracks
1828 ntr = GetNumberOfTracks();
1829 for (int itr=ntr;itr--;) {
1830 AliESDtrack *t = GetTrack(itr);
1831 if (TMath::Abs(t->Eta())>eta) continue;
1832 if (!t->IsOn(AliESDtrack::kITSin)) continue;
1833 if (useDCAFlag && t->IsOn(AliESDtrack::kMultSec)) continue;
1834 if (useV0Flag && t->IsOn(AliESDtrack::kMultInV0)) continue;
1835 if (t->IsOn(AliESDtrack::kITSpureSA)) trITSSApure++;
1841 Bool_t AliESDEvent::IsPileupFromSPDInMultBins() const {
1842 Int_t nTracklets=GetMultiplicity()->GetNumberOfTracklets();
1843 if(nTracklets<20) return IsPileupFromSPD(3,0.8);
1844 else if(nTracklets<50) return IsPileupFromSPD(4,0.8);
1845 else return IsPileupFromSPD(5,0.8);
1848 void AliESDEvent::SetTOFHeader(const AliTOFHeader *header)
1851 // Set the TOF event_time
1855 *fTOFHeader=*header;
1856 //fTOFHeader->SetName(fgkESDListName[kTOFHeader]);
1859 // for analysis of reconstructed events
1860 // when this information is not avaliable
1861 fTOFHeader = new AliTOFHeader(*header);
1862 //AddObject(fTOFHeader);
1867 AliCentrality* AliESDEvent::GetCentrality()
1869 if (!fCentrality) fCentrality = new AliCentrality();
1873 AliEventplane* AliESDEvent::GetEventplane()
1875 if (!fEventplane) fEventplane = new AliEventplane();
1879 Float_t AliESDEvent::GetVZEROEqMultiplicity(Int_t i) const
1881 // Get VZERO Multiplicity for channel i
1882 // Themethod uses the equalization factors
1883 // stored in the ESD-run object in order to
1884 // get equal multiplicities within a VZERO rins (1/8 of VZERO)
1885 if (!fESDVZERO || !fESDRun) return -1;
1888 Float_t factorSum = 0;
1889 for(Int_t j = 8*ring; j < (8*ring+8); ++j) {
1890 factorSum += fESDRun->GetVZEROEqFactors(j);
1892 Float_t factor = fESDRun->GetVZEROEqFactors(i)*8./factorSum;
1894 return (fESDVZERO->GetMultiplicity(i)/factor);