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 "AliESDTrdTrack.h"
53 #include "AliESDVertex.h"
54 #include "AliVertexerTracks.h"
55 #include "AliESDcascade.h"
56 #include "AliESDkink.h"
57 #include "AliESDtrack.h"
58 #include "AliESDHLTtrack.h"
59 #include "AliESDCaloCluster.h"
60 #include "AliESDCaloCells.h"
62 #include "AliESDFMD.h"
63 #include "AliESDVZERO.h"
64 #include "AliMultiplicity.h"
65 #include "AliRawDataErrorLog.h"
67 #include "AliESDACORDE.h"
68 #include "AliESDHLTDecision.h"
74 // here we define the names, some classes are no TNamed, therefore the classnames
76 const char* AliESDEvent::fgkESDListName[kESDListN] = {"AliESDRun",
100 "AliRawDataErrorLogs",
103 //______________________________________________________________________________
104 AliESDEvent::AliESDEvent():
106 fESDObjects(new TList()),
120 fSPDPileupVertices(0),
121 fTrkPileupVertices(0),
130 fEMCALCells(0), fPHOSCells(0),
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(esd.fESDOld ? new AliESD(*esd.fESDOld) : 0),
169 fESDFriendOld(esd.fESDFriendOld ? new AliESDfriend(*esd.fESDFriendOld) : 0),
170 fConnected(esd.fConnected),
171 fUseOwnList(esd.fUseOwnList)
174 // CKB init in the constructor list and only add here ...
179 AddObject(fESDVZERO);
180 AddObject(fESDTZERO);
181 AddObject(fTPCVertex);
182 AddObject(fSPDVertex);
183 AddObject(fPrimaryVertex);
185 AddObject(fPHOSTrigger);
186 AddObject(fEMCALTrigger);
187 AddObject(fSPDPileupVertices);
188 AddObject(fTrkPileupVertices);
190 AddObject(fMuonTracks);
191 AddObject(fPmdTracks);
192 AddObject(fTrdTracks);
194 AddObject(fCascades);
196 AddObject(fCaloClusters);
197 AddObject(fEMCALCells);
198 AddObject(fPHOSCells);
199 AddObject(fErrorLogs);
200 AddObject(fESDACORDE);
206 //______________________________________________________________________________
207 AliESDEvent & AliESDEvent::operator=(const AliESDEvent& source) {
209 // Assignment operator
211 if(&source == this) return *this;
212 AliVEvent::operator=(source);
214 // This assumes that the list is already created
215 // and that the virtual void Copy(Tobject&) function
216 // is correctly implemented in the derived class
217 // otherwise only TObject::Copy() will be used
221 if((fESDObjects->GetSize()==0)&&(source.fESDObjects->GetSize()>=kESDListN)){
222 // We cover the case that we do not yet have the
223 // standard content but the source has it
227 TIter next(source.GetList());
230 while ((its = next())) {
231 name.Form("%s", its->GetName());
232 TObject *mine = fESDObjects->FindObject(name.Data());
234 TClass* pClass=TClass::GetClass(its->ClassName());
236 AliWarning(Form("Can not find class description for entry %s (%s)\n",
237 its->ClassName(), name.Data()));
241 mine=(TObject*)pClass->New();
243 // not in this: can be added to list
244 AliWarning(Form("%s:%d Could not find %s for copying \n",
245 (char*)__FILE__,__LINE__,name.Data()));
248 if(mine->InheritsFrom("TNamed")){
249 ((TNamed*)mine)->SetName(name);
251 else if(mine->InheritsFrom("TCollection")){
252 if(mine->InheritsFrom("TClonesArray"))
253 dynamic_cast<TClonesArray*>(mine)->SetClass(dynamic_cast<TClonesArray*>(its)->GetClass());
254 dynamic_cast<TCollection*>(mine)->SetName(name);
256 AliDebug(1, Form("adding object %s of type %s", mine->GetName(), mine->ClassName()));
260 if(!its->InheritsFrom("TCollection")){
264 else if(its->InheritsFrom("TClonesArray")){
265 // Create or expand the tclonesarray pointers
266 // so we can directly copy to the object
267 TClonesArray *its_tca = (TClonesArray*)its;
268 TClonesArray *mine_tca = (TClonesArray*)mine;
270 // this leaves the capacity of the TClonesArray the same
271 // except for a factor of 2 increase when size > capacity
272 // does not release any memory occupied by the tca
273 mine_tca->ExpandCreate(its_tca->GetEntriesFast());
274 for(int i = 0;i < its_tca->GetEntriesFast();++i){
276 TObject *mine_tca_obj = mine_tca->At(i);
277 TObject *its_tca_obj = its_tca->At(i);
278 // no need to delete first
279 // pointers within the class should be handled by Copy()...
280 // Can there be Empty slots?
281 its_tca_obj->Copy(*mine_tca_obj);
285 AliWarning(Form("%s:%d cannot copy TCollection \n",
286 (char*)__FILE__,__LINE__));
290 fConnected = source.fConnected;
291 fUseOwnList = source.fUseOwnList;
298 //______________________________________________________________________________
299 AliESDEvent::~AliESDEvent()
302 // Standard destructor
305 // everthing on the list gets deleted automatically
308 if(fESDObjects&&!fConnected)
317 void AliESDEvent::Copy(TObject &obj) const {
319 // interface to TOBject::Copy
320 // Copies the content of this into obj!
321 // bascially obj = *this
323 if(this==&obj)return;
324 AliESDEvent *robj = dynamic_cast<AliESDEvent*>(&obj);
325 if(!robj)return; // not an AliESEvent
330 //______________________________________________________________________________
331 void AliESDEvent::Reset()
335 // Std content + Non std content
337 // Reset the standard contents
340 // reset for the old data without AliESDEvent...
341 if(fESDOld)fESDOld->Reset();
343 fESDFriendOld->~AliESDfriend();
344 new (fESDFriendOld) AliESDfriend();
348 if(fESDObjects->GetSize()>kESDListN){
349 // we have non std content
350 // this also covers esdfriends
351 for(int i = kESDListN;i < fESDObjects->GetSize();++i){
352 TObject *pObject = fESDObjects->At(i);
354 if(pObject->InheritsFrom(TClonesArray::Class())){
355 ((TClonesArray*)pObject)->Delete();
357 else if(!pObject->InheritsFrom(TCollection::Class())){
358 ResetWithPlacementNew(pObject);
361 AliWarning(Form("No reset for %s (%s)\n",
362 pObject->ClassName()));
369 Bool_t AliESDEvent::ResetWithPlacementNew(TObject *pObject){
370 Long_t dtoronly = TObject::GetDtorOnly();
371 TClass *pClass = TClass::GetClass(pObject->ClassName());
372 TObject::SetDtorOnly(pObject);
374 // Recreate with placement new
375 pClass->New(pObject);
376 // Restore the state.
377 TObject::SetDtorOnly((void*)dtoronly);
381 void AliESDEvent::ResetStdContent()
383 // Reset the standard contents
384 if(fESDRun) fESDRun->Reset();
385 if(fHeader) fHeader->Reset();
386 if(fESDZDC) fESDZDC->Reset();
391 // reset by callin d'to /c'tor keep the pointer
392 fESDVZERO->~AliESDVZERO();
393 new (fESDVZERO) AliESDVZERO();
396 fESDACORDE->~AliESDACORDE();
397 new (fESDACORDE) AliESDACORDE();
399 if(fESDTZERO) fESDTZERO->Reset();
400 // CKB no clear/reset implemented
402 fTPCVertex->~AliESDVertex();
403 new (fTPCVertex) AliESDVertex();
404 fTPCVertex->SetName(fgkESDListName[kTPCVertex]);
407 fSPDVertex->~AliESDVertex();
408 new (fSPDVertex) AliESDVertex();
409 fSPDVertex->SetName(fgkESDListName[kSPDVertex]);
412 fPrimaryVertex->~AliESDVertex();
413 new (fPrimaryVertex) AliESDVertex();
414 fPrimaryVertex->SetName(fgkESDListName[kPrimaryVertex]);
417 fSPDMult->~AliMultiplicity();
418 new (fSPDMult) AliMultiplicity();
420 if(fPHOSTrigger)fPHOSTrigger->Reset();
421 if(fEMCALTrigger)fEMCALTrigger->Reset();
422 if(fSPDPileupVertices)fSPDPileupVertices->Delete();
423 if(fTrkPileupVertices)fTrkPileupVertices->Delete();
424 if(fTracks)fTracks->Delete();
425 if(fMuonTracks)fMuonTracks->Delete();
426 if(fPmdTracks)fPmdTracks->Delete();
427 if(fTrdTracks)fTrdTracks->Delete();
428 if(fV0s)fV0s->Delete();
429 if(fCascades)fCascades->Delete();
430 if(fKinks)fKinks->Delete();
431 if(fCaloClusters)fCaloClusters->Delete();
432 if(fPHOSCells)fPHOSCells->DeleteContainer();
433 if(fEMCALCells)fEMCALCells->DeleteContainer();
434 if(fErrorLogs) fErrorLogs->Delete();
436 // don't reset fconnected fConnected and the list
441 Int_t AliESDEvent::AddV0(const AliESDv0 *v) {
445 TClonesArray &fv = *fV0s;
446 Int_t idx=fV0s->GetEntriesFast();
447 new(fv[idx]) AliESDv0(*v);
451 //______________________________________________________________________________
452 void AliESDEvent::Print(Option_t *) const
455 // Print header information of the event
457 printf("ESD run information\n");
458 printf("Event # in file %d Bunch crossing # %d Orbit # %d Period # %d Run # %d Trigger %lld Magnetic field %f \n",
459 GetEventNumberInFile(),
460 GetBunchCrossNumber(),
465 GetMagneticField() );
467 printf("Vertex: (%.4f +- %.4f, %.4f +- %.4f, %.4f +- %.4f) cm\n",
468 fPrimaryVertex->GetXv(), fPrimaryVertex->GetXRes(),
469 fPrimaryVertex->GetYv(), fPrimaryVertex->GetYRes(),
470 fPrimaryVertex->GetZv(), fPrimaryVertex->GetZRes());
471 printf("Mean vertex in RUN: X=%.4f Y=%.4f Z=%.4f cm\n",
472 GetDiamondX(),GetDiamondY(),GetDiamondZ());
474 printf("SPD Multiplicity. Number of tracklets %d \n",
475 fSPDMult->GetNumberOfTracklets());
476 printf("Number of pileup primary vertices reconstructed with SPD %d\n",
477 GetNumberOfPileupVerticesSPD());
478 printf("Number of pileup primary vertices reconstructed using the tracks %d\n",
479 GetNumberOfPileupVerticesTracks());
480 printf("Number of tracks: \n");
481 printf(" charged %d\n", GetNumberOfTracks());
482 printf(" muon %d\n", GetNumberOfMuonTracks());
483 printf(" pmd %d\n", GetNumberOfPmdTracks());
484 printf(" trd %d\n", GetNumberOfTrdTracks());
485 printf(" v0 %d\n", GetNumberOfV0s());
486 printf(" cascades %d\n", GetNumberOfCascades());
487 printf(" kinks %d\n", GetNumberOfKinks());
488 if(fPHOSCells)printf(" PHOSCells %d\n", fPHOSCells->GetNumberOfCells());
489 else printf(" PHOSCells not in the Event\n");
490 if(fEMCALCells)printf(" EMCALCells %d\n", fEMCALCells->GetNumberOfCells());
491 else printf(" EMCALCells not in the Event\n");
492 printf(" CaloClusters %d\n", GetNumberOfCaloClusters());
493 printf(" FMD %s\n", (fESDFMD ? "yes" : "no"));
494 printf(" VZERO %s\n", (fESDVZERO ? "yes" : "no"));
495 TObject* pHLTDecision=GetHLTTriggerDecision();
496 printf("HLT trigger decision: %s\n", pHLTDecision?pHLTDecision->GetOption():"not available");
497 if (pHLTDecision) pHLTDecision->Print("compact");
502 void AliESDEvent::SetESDfriend(const AliESDfriend *ev) const {
504 // Attaches the complementary info to the ESD
508 // to be sure that we set the tracks also
509 // in case of old esds
510 // if(fESDOld)CopyFromOldESD();
512 Int_t ntrk=ev->GetNumberOfTracks();
514 for (Int_t i=0; i<ntrk; i++) {
515 const AliESDfriendTrack *f=ev->GetTrack(i);
516 GetTrack(i)->SetFriendTrack(f);
520 Bool_t AliESDEvent::RemoveKink(Int_t rm) const {
521 // ---------------------------------------------------------
522 // Remove a kink candidate and references to it from ESD,
523 // if this candidate does not come from a reconstructed decay
524 // Not yet implemented...
525 // ---------------------------------------------------------
526 Int_t last=GetNumberOfKinks()-1;
527 if ((rm<0)||(rm>last)) return kFALSE;
532 Bool_t AliESDEvent::RemoveV0(Int_t rm) const {
533 // ---------------------------------------------------------
534 // Remove a V0 candidate and references to it from ESD,
535 // if this candidate does not come from a reconstructed decay
536 // ---------------------------------------------------------
537 Int_t last=GetNumberOfV0s()-1;
538 if ((rm<0)||(rm>last)) return kFALSE;
540 AliESDv0 *v0=GetV0(rm);
541 Int_t idxP=v0->GetPindex(), idxN=v0->GetNindex();
544 Int_t lastIdxP=v0->GetPindex(), lastIdxN=v0->GetNindex();
548 // Check if this V0 comes from a reconstructed decay
549 Int_t ncs=GetNumberOfCascades();
550 for (Int_t n=0; n<ncs; n++) {
551 AliESDcascade *cs=GetCascade(n);
553 Int_t csIdxP=cs->GetPindex();
554 Int_t csIdxN=cs->GetNindex();
557 if (idxN==csIdxN) return kFALSE;
559 if (csIdxP==lastIdxP)
560 if (csIdxN==lastIdxN) used++;
563 //Replace the removed V0 with the last V0
564 TClonesArray &a=*fV0s;
565 delete a.RemoveAt(rm);
567 if (rm==last) return kTRUE;
569 //v0 is pointing to the last V0 candidate...
570 new (a[rm]) AliESDv0(*v0);
571 delete a.RemoveAt(last);
573 if (!used) return kTRUE;
576 // Remap the indices of the daughters of reconstructed decays
577 for (Int_t n=0; n<ncs; n++) {
578 AliESDcascade *cs=GetCascade(n);
581 Int_t csIdxP=cs->GetPindex();
582 Int_t csIdxN=cs->GetNindex();
584 if (csIdxP==lastIdxP)
585 if (csIdxN==lastIdxN) {
586 cs->AliESDv0::SetIndex(1,idxP);
587 cs->AliESDv0::SetIndex(0,idxN);
589 if (!used) return kTRUE;
596 Bool_t AliESDEvent::RemoveTrack(Int_t rm) const {
597 // ---------------------------------------------------------
598 // Remove a track and references to it from ESD,
599 // if this track does not come from a reconstructed decay
600 // ---------------------------------------------------------
601 Int_t last=GetNumberOfTracks()-1;
602 if ((rm<0)||(rm>last)) return kFALSE;
606 // Check if this track comes from the reconstructed primary vertices
607 if (fTPCVertex && fTPCVertex->GetStatus()) {
608 UShort_t *primIdx=fTPCVertex->GetIndices();
609 Int_t n=fTPCVertex->GetNIndices();
611 Int_t idx=Int_t(primIdx[n]);
612 if (rm==idx) return kFALSE;
613 if (idx==last) used++;
616 if (fPrimaryVertex && fPrimaryVertex->GetStatus()) {
617 UShort_t *primIdx=fPrimaryVertex->GetIndices();
618 Int_t n=fPrimaryVertex->GetNIndices();
620 Int_t idx=Int_t(primIdx[n]);
621 if (rm==idx) return kFALSE;
622 if (idx==last) used++;
626 // Check if this track comes from a reconstructed decay
627 Int_t nv0=GetNumberOfV0s();
628 for (Int_t n=0; n<nv0; n++) {
629 AliESDv0 *v0=GetV0(n);
631 Int_t idx=v0->GetNindex();
632 if (rm==idx) return kFALSE;
633 if (idx==last) used++;
636 if (rm==idx) return kFALSE;
637 if (idx==last) used++;
640 Int_t ncs=GetNumberOfCascades();
641 for (Int_t n=0; n<ncs; n++) {
642 AliESDcascade *cs=GetCascade(n);
644 Int_t idx=cs->GetIndex();
645 if (rm==idx) return kFALSE;
646 if (idx==last) used++;
650 if (rm==idx) return kFALSE;
651 if (idx==last) used++;
654 if (rm==idx) return kFALSE;
655 if (idx==last) used++;
658 Int_t nkn=GetNumberOfKinks();
659 for (Int_t n=0; n<nkn; n++) {
660 AliESDkink *kn=GetKink(n);
662 Int_t idx=kn->GetIndex(0);
663 if (rm==idx) return kFALSE;
664 if (idx==last) used++;
667 if (rm==idx) return kFALSE;
668 if (idx==last) used++;
671 // Check if this track is associated with a CaloCluster
672 Int_t ncl=GetNumberOfCaloClusters();
673 for (Int_t n=0; n<ncl; n++) {
674 AliESDCaloCluster *cluster=GetCaloCluster(n);
675 TArrayI *arr=cluster->GetTracksMatched();
676 Int_t s=arr->GetSize();
678 Int_t idx=arr->At(s);
679 if (rm==idx) return kFALSE;
680 if (idx==last) used++;
686 //Replace the removed track with the last track
687 TClonesArray &a=*fTracks;
688 delete a.RemoveAt(rm);
690 if (rm==last) return kTRUE;
692 AliESDtrack *t=GetTrack(last);
694 new (a[rm]) AliESDtrack(*t);
695 delete a.RemoveAt(last);
698 if (!used) return kTRUE;
701 // Remap the indices of the tracks used for the primary vertex reconstruction
702 if (fTPCVertex && fTPCVertex->GetStatus()) {
703 UShort_t *primIdx=fTPCVertex->GetIndices();
704 Int_t n=fTPCVertex->GetNIndices();
706 Int_t idx=Int_t(primIdx[n]);
708 primIdx[n]=Short_t(rm);
710 if (!used) return kTRUE;
714 if (fPrimaryVertex && fPrimaryVertex->GetStatus()) {
715 UShort_t *primIdx=fPrimaryVertex->GetIndices();
716 Int_t n=fPrimaryVertex->GetNIndices();
718 Int_t idx=Int_t(primIdx[n]);
720 primIdx[n]=Short_t(rm);
722 if (!used) return kTRUE;
727 // Remap the indices of the daughters of reconstructed decays
728 for (Int_t n=0; n<nv0; n++) {
729 AliESDv0 *v0=GetV0(n);
730 if (v0->GetIndex(0)==last) {
733 if (!used) return kTRUE;
735 if (v0->GetIndex(1)==last) {
738 if (!used) return kTRUE;
742 for (Int_t n=0; n<ncs; n++) {
743 AliESDcascade *cs=GetCascade(n);
744 if (cs->GetIndex()==last) {
747 if (!used) return kTRUE;
750 if (v0->GetIndex(0)==last) {
753 if (!used) return kTRUE;
755 if (v0->GetIndex(1)==last) {
758 if (!used) return kTRUE;
762 for (Int_t n=0; n<nkn; n++) {
763 AliESDkink *kn=GetKink(n);
764 if (kn->GetIndex(0)==last) {
767 if (!used) return kTRUE;
769 if (kn->GetIndex(1)==last) {
772 if (!used) return kTRUE;
776 // Remap the indices of the tracks accosicated with CaloClusters
777 for (Int_t n=0; n<ncl; n++) {
778 AliESDCaloCluster *cluster=GetCaloCluster(n);
779 TArrayI *arr=cluster->GetTracksMatched();
780 Int_t s=arr->GetSize();
782 Int_t idx=arr->At(s);
786 if (!used) return kTRUE;
795 Bool_t AliESDEvent::Clean(Float_t *cleanPars) {
797 // Remove the data which are not needed for the physics analysis.
799 // 1) Cleaning the V0 candidates
800 // ---------------------------
801 // If the cosine of the V0 pointing angle "csp" and
802 // the DCA between the daughter tracks "dca" does not satisfy
805 // csp > cleanPars[1] + dca/cleanPars[0]*(1.- cleanPars[1])
807 // an attempt to remove this V0 candidate from ESD is made.
809 // The V0 candidate gets removed if it does not belong to any
810 // recosntructed cascade decay
812 // 12.11.2007, optimal values: cleanPars[0]=0.5, cleanPars[1]=0.999
814 // 2) Cleaning the tracks
815 // ----------------------
816 // If track's transverse parameter is larger than cleanPars[2]
818 // track's longitudinal parameter is larger than cleanPars[3]
819 // an attempt to remove this track from ESD is made.
821 // The track gets removed if it does not come
822 // from a reconstructed decay
826 Float_t dcaMax=cleanPars[0];
827 Float_t cspMin=cleanPars[1];
829 Int_t nV0s=GetNumberOfV0s();
830 for (Int_t i=nV0s-1; i>=0; i--) {
831 AliESDv0 *v0=GetV0(i);
833 Float_t dca=v0->GetDcaV0Daughters();
834 Float_t csp=v0->GetV0CosineOfPointingAngle();
835 Float_t cspcut=cspMin + dca/dcaMax*(1.-cspMin);
836 if (csp > cspcut) continue;
837 if (RemoveV0(i)) rc=kTRUE;
841 Float_t dmax=cleanPars[2], zmax=cleanPars[3];
843 const AliESDVertex *vertex=GetPrimaryVertexSPD();
844 Bool_t vtxOK=vertex->GetStatus();
846 Int_t nTracks=GetNumberOfTracks();
847 for (Int_t i=nTracks-1; i>=0; i--) {
848 AliESDtrack *track=GetTrack(i);
849 Float_t xy,z; track->GetImpactParameters(xy,z);
850 if ((TMath::Abs(xy) > dmax) || (vtxOK && (TMath::Abs(z) > zmax))) {
851 if (RemoveTrack(i)) rc=kTRUE;
858 Char_t AliESDEvent::AddPileupVertexSPD(const AliESDVertex *vtx)
860 // Add a pileup primary vertex reconstructed with SPD
861 TClonesArray &ftr = *fSPDPileupVertices;
862 Char_t n=Char_t(ftr.GetEntriesFast());
863 AliESDVertex *vertex = new(ftr[n]) AliESDVertex(*vtx);
868 Char_t AliESDEvent::AddPileupVertexTracks(const AliESDVertex *vtx)
870 // Add a pileup primary vertex reconstructed with SPD
871 TClonesArray &ftr = *fTrkPileupVertices;
872 Char_t n=Char_t(ftr.GetEntriesFast());
873 AliESDVertex *vertex = new(ftr[n]) AliESDVertex(*vtx);
878 Int_t AliESDEvent::AddTrack(const AliESDtrack *t)
881 TClonesArray &ftr = *fTracks;
882 AliESDtrack * track = new(ftr[fTracks->GetEntriesFast()])AliESDtrack(*t);
883 track->SetID(fTracks->GetEntriesFast()-1);
884 return track->GetID();
887 void AliESDEvent::AddMuonTrack(const AliESDMuonTrack *t)
889 TClonesArray &fmu = *fMuonTracks;
890 new(fmu[fMuonTracks->GetEntriesFast()]) AliESDMuonTrack(*t);
893 void AliESDEvent::AddPmdTrack(const AliESDPmdTrack *t)
895 TClonesArray &fpmd = *fPmdTracks;
896 new(fpmd[fPmdTracks->GetEntriesFast()]) AliESDPmdTrack(*t);
899 void AliESDEvent::AddTrdTrack(const AliESDTrdTrack *t)
901 TClonesArray &ftrd = *fTrdTracks;
902 new(ftrd[fTrdTracks->GetEntriesFast()]) AliESDTrdTrack(*t);
908 Int_t AliESDEvent::AddKink(const AliESDkink *c)
911 TClonesArray &fk = *fKinks;
912 AliESDkink * kink = new(fk[fKinks->GetEntriesFast()]) AliESDkink(*c);
913 kink->SetID(fKinks->GetEntriesFast()); // CKB different from the other imps..
914 return fKinks->GetEntriesFast()-1;
918 void AliESDEvent::AddCascade(const AliESDcascade *c)
920 TClonesArray &fc = *fCascades;
921 new(fc[fCascades->GetEntriesFast()]) AliESDcascade(*c);
925 Int_t AliESDEvent::AddCaloCluster(const AliESDCaloCluster *c)
928 TClonesArray &fc = *fCaloClusters;
929 AliESDCaloCluster *clus = new(fc[fCaloClusters->GetEntriesFast()]) AliESDCaloCluster(*c);
930 clus->SetID(fCaloClusters->GetEntriesFast()-1);
931 return fCaloClusters->GetEntriesFast()-1;
935 void AliESDEvent::AddRawDataErrorLog(const AliRawDataErrorLog *log) const {
936 TClonesArray &errlogs = *fErrorLogs;
937 new(errlogs[errlogs.GetEntriesFast()]) AliRawDataErrorLog(*log);
940 void AliESDEvent::SetPrimaryVertexTPC(const AliESDVertex *vertex)
942 // Set the TPC vertex
943 // use already allocated space
945 *fTPCVertex = *vertex;
946 fTPCVertex->SetName(fgkESDListName[kTPCVertex]);
950 void AliESDEvent::SetPrimaryVertexSPD(const AliESDVertex *vertex)
952 // Set the SPD vertex
953 // use already allocated space
955 *fSPDVertex = *vertex;
956 fSPDVertex->SetName(fgkESDListName[kSPDVertex]);
960 void AliESDEvent::SetPrimaryVertexTracks(const AliESDVertex *vertex)
962 // Set the primary vertex reconstructed using he ESD tracks.
963 // use already allocated space
965 *fPrimaryVertex = *vertex;
966 fPrimaryVertex->SetName(fgkESDListName[kPrimaryVertex]);
970 const AliESDVertex * AliESDEvent::GetPrimaryVertex() const
973 // Get the "best" available reconstructed primary vertex.
976 if (fPrimaryVertex->GetStatus()) return fPrimaryVertex;
979 if (fSPDVertex->GetStatus()) return fSPDVertex;
981 if(fTPCVertex) return fTPCVertex;
983 AliWarning("No primary vertex available. Returning the \"default\"...");
987 AliESDVertex * AliESDEvent::PrimaryVertexTracksUnconstrained() const
990 // Removes diamond constraint from fPrimaryVertex (reconstructed with tracks)
991 // Returns a AliESDVertex which has to be deleted by the user
993 if(!fPrimaryVertex) {
994 AliWarning("No primary vertex from tracks available.");
997 if(!fPrimaryVertex->GetStatus()) {
998 AliWarning("No primary vertex from tracks available.");
1002 AliVertexerTracks vertexer(GetMagneticField());
1003 Float_t diamondxyz[3]={(Float_t)GetDiamondX(),(Float_t)GetDiamondY(),0.};
1004 Float_t diamondcovxy[3]; GetDiamondCovXY(diamondcovxy);
1005 Float_t diamondcov[6]={diamondcovxy[0],diamondcovxy[1],diamondcovxy[2],0.,0.,7.};
1006 AliESDVertex *vertex =
1007 (AliESDVertex*)vertexer.RemoveConstraintFromVertex(fPrimaryVertex,diamondxyz,diamondcov);
1012 void AliESDEvent::SetMultiplicity(const AliMultiplicity *mul)
1014 // Set the SPD Multiplicity
1021 void AliESDEvent::SetFMDData(AliESDFMD * obj)
1023 // use already allocated space
1029 void AliESDEvent::SetVZEROData(AliESDVZERO * obj)
1031 // use already allocated space
1036 void AliESDEvent::SetACORDEData(AliESDACORDE * obj)
1043 void AliESDEvent::GetESDfriend(AliESDfriend *ev) const
1046 // Extracts the complementary info from the ESD
1050 Int_t ntrk=GetNumberOfTracks();
1052 for (Int_t i=0; i<ntrk; i++) {
1053 AliESDtrack *t=GetTrack(i);
1054 const AliESDfriendTrack *f=t->GetFriendTrack();
1057 t->ReleaseESDfriendTrack();// Not to have two copies of "friendTrack"
1061 AliESDfriend *fr = (AliESDfriend*)(const_cast<AliESDEvent*>(this)->FindListObject("AliESDfriend"));
1062 if (fr) ev->SetVZEROfriend(fr->GetVZEROfriend());
1065 void AliESDEvent::AddObject(TObject* obj)
1067 // Add an object to the list of object.
1068 // Please be aware that in order to increase performance you should
1069 // refrain from using TObjArrays (if possible). Use TClonesArrays, instead.
1070 fESDObjects->SetOwner(kTRUE);
1071 fESDObjects->AddLast(obj);
1075 void AliESDEvent::GetStdContent()
1077 // set pointers for standard content
1078 // get by name much safer and not a big overhead since not called very often
1080 fESDRun = (AliESDRun*)fESDObjects->FindObject(fgkESDListName[kESDRun]);
1081 fHeader = (AliESDHeader*)fESDObjects->FindObject(fgkESDListName[kHeader]);
1082 fESDZDC = (AliESDZDC*)fESDObjects->FindObject(fgkESDListName[kESDZDC]);
1083 fESDFMD = (AliESDFMD*)fESDObjects->FindObject(fgkESDListName[kESDFMD]);
1084 fESDVZERO = (AliESDVZERO*)fESDObjects->FindObject(fgkESDListName[kESDVZERO]);
1085 fESDTZERO = (AliESDTZERO*)fESDObjects->FindObject(fgkESDListName[kESDTZERO]);
1086 fTPCVertex = (AliESDVertex*)fESDObjects->FindObject(fgkESDListName[kTPCVertex]);
1087 fSPDVertex = (AliESDVertex*)fESDObjects->FindObject(fgkESDListName[kSPDVertex]);
1088 fPrimaryVertex = (AliESDVertex*)fESDObjects->FindObject(fgkESDListName[kPrimaryVertex]);
1089 fSPDMult = (AliMultiplicity*)fESDObjects->FindObject(fgkESDListName[kSPDMult]);
1090 fPHOSTrigger = (AliESDCaloTrigger*)fESDObjects->FindObject(fgkESDListName[kPHOSTrigger]);
1091 fEMCALTrigger = (AliESDCaloTrigger*)fESDObjects->FindObject(fgkESDListName[kEMCALTrigger]);
1092 fSPDPileupVertices = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kSPDPileupVertices]);
1093 fTrkPileupVertices = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kTrkPileupVertices]);
1094 fTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kTracks]);
1095 fMuonTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kMuonTracks]);
1096 fPmdTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kPmdTracks]);
1097 fTrdTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kTrdTracks]);
1098 fV0s = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kV0s]);
1099 fCascades = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kCascades]);
1100 fKinks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kKinks]);
1101 fCaloClusters = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kCaloClusters]);
1102 fEMCALCells = (AliESDCaloCells*)fESDObjects->FindObject(fgkESDListName[kEMCALCells]);
1103 fPHOSCells = (AliESDCaloCells*)fESDObjects->FindObject(fgkESDListName[kPHOSCells]);
1104 fErrorLogs = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kErrorLogs]);
1105 fESDACORDE = (AliESDACORDE*)fESDObjects->FindObject(fgkESDListName[kESDACORDE]);
1109 void AliESDEvent::SetStdNames(){
1110 // Set the names of the standard contents
1112 if(fESDObjects->GetEntries()>=kESDListN){
1113 for(int i = 0;i < fESDObjects->GetEntries() && i<kESDListN;i++){
1114 TObject *fObj = fESDObjects->At(i);
1115 if(fObj->InheritsFrom("TNamed")){
1116 ((TNamed*)fObj)->SetName(fgkESDListName[i]);
1118 else if(fObj->InheritsFrom("TClonesArray")){
1119 ((TClonesArray*)fObj)->SetName(fgkESDListName[i]);
1124 AliWarning("Std Entries missing");
1129 void AliESDEvent::CreateStdContent(Bool_t bUseThisList){
1130 fUseOwnList = bUseThisList;
1134 void AliESDEvent::CreateStdContent()
1136 // create the standard AOD content and set pointers
1138 // create standard objects and add them to the TList of objects
1139 AddObject(new AliESDRun());
1140 AddObject(new AliESDHeader());
1141 AddObject(new AliESDZDC());
1142 AddObject(new AliESDFMD());
1143 AddObject(new AliESDVZERO());
1144 AddObject(new AliESDTZERO());
1145 AddObject(new AliESDVertex());
1146 AddObject(new AliESDVertex());
1147 AddObject(new AliESDVertex());
1148 AddObject(new AliMultiplicity());
1149 AddObject(new AliESDCaloTrigger());
1150 AddObject(new AliESDCaloTrigger());
1151 AddObject(new TClonesArray("AliESDVertex",0));
1152 AddObject(new TClonesArray("AliESDVertex",0));
1153 AddObject(new TClonesArray("AliESDtrack",0));
1154 AddObject(new TClonesArray("AliESDMuonTrack",0));
1155 AddObject(new TClonesArray("AliESDPmdTrack",0));
1156 AddObject(new TClonesArray("AliESDTrdTrack",0));
1157 AddObject(new TClonesArray("AliESDv0",0));
1158 AddObject(new TClonesArray("AliESDcascade",0));
1159 AddObject(new TClonesArray("AliESDkink",0));
1160 AddObject(new TClonesArray("AliESDCaloCluster",0));
1161 AddObject(new AliESDCaloCells());
1162 AddObject(new AliESDCaloCells());
1163 AddObject(new TClonesArray("AliRawDataErrorLog",0));
1164 AddObject(new AliESDACORDE());
1166 // check the order of the indices against enum...
1170 // read back pointers
1174 TObject* AliESDEvent::FindListObject(const char *name) const {
1176 // Find object with name "name" in the list of branches
1179 return fESDObjects->FindObject(name);
1184 Int_t AliESDEvent::GetPHOSClusters(TRefArray *clusters) const
1186 // fills the provided TRefArray with all found phos clusters
1190 AliESDCaloCluster *cl = 0;
1191 for (Int_t i = 0; i < GetNumberOfCaloClusters(); i++) {
1193 if ( (cl = GetCaloCluster(i)) ) {
1196 AliDebug(1,Form("IsPHOS cluster %d Size: %d \n",i,clusters->GetEntriesFast()));
1200 return clusters->GetEntriesFast();
1203 Int_t AliESDEvent::GetEMCALClusters(TRefArray *clusters) const
1205 // fills the provided TRefArray with all found emcal clusters
1209 AliESDCaloCluster *cl = 0;
1210 for (Int_t i = 0; i < GetNumberOfCaloClusters(); i++) {
1212 if ( (cl = GetCaloCluster(i)) ) {
1215 AliDebug(1,Form("IsEMCAL cluster %d Size: %d \n",i,clusters->GetEntriesFast()));
1219 return clusters->GetEntriesFast();
1222 void AliESDEvent::WriteToTree(TTree* tree) const {
1223 // Book the branches as in TTree::Branch(TCollection*)
1224 // but add a "." at the end of top level branches which are
1225 // not a TClonesArray
1229 TIter next(fESDObjects);
1230 const Int_t kSplitlevel = 99; // default value in TTree::Branch()
1231 const Int_t kBufsize = 32000; // default value in TTree::Branch()
1234 while ((obj = next())) {
1235 branchname.Form("%s", obj->GetName());
1236 if(branchname.CompareTo("AliESDfriend")==0)branchname = "ESDfriend.";
1237 if ((kSplitlevel > 1) && !obj->InheritsFrom(TClonesArray::Class())) {
1238 if(!branchname.EndsWith("."))branchname += ".";
1240 if (!tree->FindBranch(branchname)) {
1241 // For the custom streamer to be called splitlevel
1242 // has to be negative, only needed for HLT
1243 Int_t splitLevel = (TString(obj->ClassName()) == "AliHLTGlobalTriggerDecision") ? -1 : kSplitlevel - 1;
1244 tree->Bronch(branchname, obj->ClassName(), fESDObjects->GetObjectRef(obj),kBufsize, splitLevel);
1250 void AliESDEvent::ReadFromTree(TTree *tree, Option_t* opt){
1252 // Connect the ESDEvent to a tree
1255 AliWarning("AliESDEvent::ReadFromTree() Zero Pointer to Tree \n");
1259 if(!tree->GetTree())tree->LoadTree(0);
1261 // if we find the "ESD" branch on the tree we do have the old structure
1262 if(tree->GetBranch("ESD")) {
1263 char ** address = (char **)(tree->GetBranch("ESD")->GetAddress());
1264 // do we have the friend branch
1265 TBranch * esdFB = tree->GetBranch("ESDfriend.");
1266 char ** addressF = 0;
1267 if(esdFB)addressF = (char **)(esdFB->GetAddress());
1269 AliInfo("AliESDEvent::ReadFromTree() Reading old Tree");
1270 tree->SetBranchAddress("ESD", &fESDOld);
1272 tree->SetBranchAddress("ESDfriend.",&fESDFriendOld);
1275 AliInfo("AliESDEvent::ReadFromTree() Reading old Tree");
1276 AliInfo("Branch already connected. Using existing branch address.");
1277 fESDOld = (AliESD*) (*address);
1278 // addressF can still be 0, since branch needs to switched on
1279 if(addressF)fESDFriendOld = (AliESDfriend*) (*addressF);
1282 // have already connected the old ESD structure... ?
1283 // reuse also the pointer of the AlliESDEvent
1284 // otherwise create new ones
1285 TList* connectedList = (TList*) (tree->GetUserInfo()->FindObject("ESDObjectsConnectedToTree"));
1288 // If connected use the connected list of objects
1289 if(fESDObjects!= connectedList){
1290 // protect when called twice
1291 fESDObjects->Delete();
1292 fESDObjects = connectedList;
1297 // The pointer to the friend changes when called twice via InitIO
1298 // since AliESDEvent is deleted
1299 TObject* oldf = FindListObject("AliESDfriend");
1302 newf = (TObject*)*addressF;
1304 if(newf!=0&&oldf!=newf){
1305 // remove the old reference
1306 // Should we also delete it? Or is this handled in TTree I/O
1307 // since it is created by the first SetBranchAddress
1308 fESDObjects->Remove(oldf);
1310 fESDObjects->Add(newf);
1317 CreateStdContent(); // create for copy
1318 // if we have the esdfriend add it, so we always can access it via the userinfo
1319 if(fESDFriendOld)AddObject(fESDFriendOld);
1320 // we are not owner of the list objects
1321 // must not delete it
1322 fESDObjects->SetOwner(kTRUE);
1323 fESDObjects->SetName("ESDObjectsConnectedToTree");
1324 tree->GetUserInfo()->Add(fESDObjects);
1332 // Try to find AliESDEvent
1333 AliESDEvent *esdEvent = 0;
1334 esdEvent = (AliESDEvent*)tree->GetTree()->GetUserInfo()->FindObject("AliESDEvent");
1336 // Check if already connected to tree
1338 TList* connectedList = (TList*) (tree->GetUserInfo()->FindObject("ESDObjectsConnectedToTree"));
1341 if (connectedList && (strcmp(opt, "reconnect"))) {
1342 // If connected use the connected list if objects
1343 fESDObjects->Delete();
1344 fESDObjects = connectedList;
1351 // prevent a memory leak when reading back the TList
1352 // if (!(strcmp(opt, "reconnect"))) fESDObjects->Delete();
1355 // create a new TList from the UserInfo TList...
1356 // copy constructor does not work...
1357 fESDObjects = (TList*)(esdEvent->GetList()->Clone());
1358 fESDObjects->SetOwner(kTRUE);
1360 else if ( fESDObjects->GetEntries()==0){
1361 // at least create the std content if we want to read to our list
1366 // we only need new things in the list if we do no already have it..
1367 // TODO just add new entries
1369 if(fESDObjects->GetEntries()<kESDListN){
1370 AliWarning(Form("AliESDEvent::ReadFromTree() TList contains less than the standard contents %d < %d \n",
1371 fESDObjects->GetEntries(),kESDListN));
1373 // set the branch addresses
1374 TIter next(fESDObjects);
1376 while((el=(TNamed*)next())){
1377 TString bname(el->GetName());
1378 if(bname.CompareTo("AliESDfriend")==0)
1380 // AliESDfriend does not have a name ...
1381 tree->SetBranchAddress("ESDfriend.",fESDObjects->GetObjectRef(el));
1384 // check if branch exists under this Name
1385 TBranch *br = tree->GetBranch(bname.Data());
1387 tree->SetBranchAddress(bname.Data(),fESDObjects->GetObjectRef(el));
1390 br = tree->GetBranch(Form("%s.",bname.Data()));
1392 tree->SetBranchAddress(Form("%s.",bname.Data()),fESDObjects->GetObjectRef(el));
1395 AliWarning(Form("AliESDEvent::ReadFromTree() No Branch found with Name %s or %s.",bname.Data(),bname.Data()));
1402 // when reading back we are not owner of the list
1403 // must not delete it
1404 fESDObjects->SetOwner(kTRUE);
1405 fESDObjects->SetName("ESDObjectsConnectedToTree");
1406 // we are not owner of the list objects
1407 // must not delete it
1408 tree->GetUserInfo()->Add(fESDObjects);
1409 tree->GetUserInfo()->SetOwner(kFALSE);
1413 // we can't get the list from the user data, create standard content
1414 // and set it by hand (no ESDfriend at the moment
1416 TIter next(fESDObjects);
1418 while((el=(TNamed*)next())){
1419 TString bname(el->GetName());
1420 TBranch *br = tree->GetBranch(bname.Data());
1422 tree->SetBranchAddress(bname.Data(),fESDObjects->GetObjectRef(el));
1425 br = tree->GetBranch(Form("%s.",bname.Data()));
1427 tree->SetBranchAddress(Form("%s.",bname.Data()),fESDObjects->GetObjectRef(el));
1432 // when reading back we are not owner of the list
1433 // must not delete it
1434 fESDObjects->SetOwner(kTRUE);
1439 void AliESDEvent::CopyFromOldESD()
1441 // Method which copies over everthing from the old esd structure to the
1446 SetRunNumber(fESDOld->GetRunNumber());
1447 SetPeriodNumber(fESDOld->GetPeriodNumber());
1448 SetMagneticField(fESDOld->GetMagneticField());
1450 // leave out diamond ...
1451 // SetDiamond(const AliESDVertex *vertex) { fESDRun->SetDiamond(vertex);}
1454 SetTriggerMask(fESDOld->GetTriggerMask());
1455 SetOrbitNumber(fESDOld->GetOrbitNumber());
1456 SetTimeStamp(fESDOld->GetTimeStamp());
1457 SetEventType(fESDOld->GetEventType());
1458 SetEventNumberInFile(fESDOld->GetEventNumberInFile());
1459 SetBunchCrossNumber(fESDOld->GetBunchCrossNumber());
1460 SetTriggerCluster(fESDOld->GetTriggerCluster());
1464 SetZDC(fESDOld->GetZDCN1Energy(),
1465 fESDOld->GetZDCP1Energy(),
1466 fESDOld->GetZDCEMEnergy(),
1468 fESDOld->GetZDCN2Energy(),
1469 fESDOld->GetZDCP2Energy(),
1470 fESDOld->GetZDCParticipants(),
1480 if(fESDOld->GetFMDData())SetFMDData(fESDOld->GetFMDData());
1484 SetT0zVertex(fESDOld->GetT0zVertex());
1485 SetT0(fESDOld->GetT0());
1489 if (fESDOld->GetVZEROData()) SetVZEROData(fESDOld->GetVZEROData());
1491 if(fESDOld->GetVertex())SetPrimaryVertexSPD(fESDOld->GetVertex());
1493 if(fESDOld->GetPrimaryVertex())SetPrimaryVertexTracks(fESDOld->GetPrimaryVertex());
1495 if(fESDOld->GetMultiplicity())SetMultiplicity(fESDOld->GetMultiplicity());
1497 for(int i = 0;i<fESDOld->GetNumberOfTracks();i++){
1498 AddTrack(fESDOld->GetTrack(i));
1501 for(int i = 0;i<fESDOld->GetNumberOfMuonTracks();i++){
1502 AddMuonTrack(fESDOld->GetMuonTrack(i));
1505 for(int i = 0;i<fESDOld->GetNumberOfPmdTracks();i++){
1506 AddPmdTrack(fESDOld->GetPmdTrack(i));
1509 for(int i = 0;i<fESDOld->GetNumberOfTrdTracks();i++){
1510 AddTrdTrack(fESDOld->GetTrdTrack(i));
1513 for(int i = 0;i<fESDOld->GetNumberOfV0s();i++){
1514 AddV0(fESDOld->GetV0(i));
1517 for(int i = 0;i<fESDOld->GetNumberOfCascades();i++){
1518 AddCascade(fESDOld->GetCascade(i));
1521 for(int i = 0;i<fESDOld->GetNumberOfKinks();i++){
1522 AddKink(fESDOld->GetKink(i));
1526 for(int i = 0;i<fESDOld->GetNumberOfCaloClusters();i++){
1527 AddCaloCluster(fESDOld->GetCaloCluster(i));
1533 Bool_t AliESDEvent::IsEventSelected(const char *trigExpr) const
1535 // Check if the event satisfies the trigger
1536 // selection expression trigExpr.
1537 // trigExpr can be any logical expression
1538 // of the trigger classes defined in AliESDRun
1539 // In case of wrong syntax return kTRUE.
1541 TString expr(trigExpr);
1542 if (expr.IsNull()) return kTRUE;
1544 ULong64_t mask = GetTriggerMask();
1545 for(Int_t itrig = 0; itrig < AliESDRun::kNTriggerClasses; itrig++) {
1546 if (mask & (1ull << itrig)) {
1547 expr.ReplaceAll(GetESDRun()->GetTriggerClass(itrig),"1");
1550 expr.ReplaceAll(GetESDRun()->GetTriggerClass(itrig),"0");
1555 if ((gROOT->ProcessLineFast(expr.Data(),&error) == 0) &&
1556 (error == TInterpreter::kNoError)) {
1564 TObject* AliESDEvent::GetHLTTriggerDecision() const
1566 // get the HLT trigger decission object
1568 // cast away const'nes because the FindListObject method
1570 AliESDEvent* pNonConst=const_cast<AliESDEvent*>(this);
1571 return pNonConst->FindListObject("HLTGlobalTrigger");
1574 TString AliESDEvent::GetHLTTriggerDescription() const
1576 // get the HLT trigger decission description
1577 TString description;
1578 TObject* pDecision=GetHLTTriggerDecision();
1580 description=pDecision->GetTitle();
1586 Bool_t AliESDEvent::IsHLTTriggerFired(const char* name) const
1588 // get the HLT trigger decission description
1589 TObject* pDecision=GetHLTTriggerDecision();
1590 if (!pDecision) return kFALSE;
1592 Option_t* option=pDecision->GetOption();
1593 if (option==NULL || *option!='1') return kFALSE;
1596 TString description=GetHLTTriggerDescription();
1597 Int_t index=description.Index(name);
1598 if (index<0) return kFALSE;
1599 index+=strlen(name);
1600 if (index>=description.Length()) return kFALSE;
1601 if (description[index]!=0 && description[index]!=' ') return kFALSE;
1606 Bool_t AliESDEvent::IsPileupFromSPD(Int_t ncont, Double_t distz, Double_t nSigmaDeltaZ, Double_t nSigmaXY, Int_t option) const {
1608 // This function checks if there was a pile up
1609 // reconstructed with SPD
1611 Double_t diamx= GetDiamondX();
1612 Double_t diamsigma2x= GetSigma2DiamondX();
1613 Double_t diamy= GetDiamondY();
1614 Double_t diamsigma2y= GetSigma2DiamondY();
1616 Double_t sigmax= TMath::Sqrt(diamsigma2x);
1617 Double_t sigmay= TMath::Sqrt(diamsigma2y);
1619 Double_t z1=fSPDVertex->GetZ();
1620 Int_t nc1=fSPDVertex->GetNContributors();
1621 if(nc1<1) return kFALSE;
1622 Int_t nPileVert=GetNumberOfPileupVerticesSPD();
1623 if(nPileVert==0) return kFALSE;
1624 for(Int_t i=0; i<nPileVert;i++){
1625 const AliESDVertex* pv=GetPileupVertexSPD(i);
1626 Double_t z2=pv->GetZ();
1627 Double_t x2=pv->GetX();
1628 Double_t y2=pv->GetY();
1629 Int_t nc2=pv->GetNContributors();
1630 Double_t distanceZ=TMath::Abs(z2-z1);
1631 Double_t distanceX=TMath::Abs(x2-diamx);
1632 Double_t distanceY=TMath::Abs(y2-diamy);
1633 Double_t errzDist=0.;
1634 Double_t errxDist=0.;
1635 Double_t erryDist=0.;
1637 Double_t ez1=fSPDVertex->GetZRes();
1638 Double_t ez2=pv->GetZRes();
1639 errzDist=TMath::Sqrt(ez1*ez1+ez2*ez2);
1641 Double_t resol1=-75.6+834.6/TMath::Sqrt((Double_t)nc1);
1643 Double_t resol2=-75.6+834.6/TMath::Sqrt((Double_t)nc2);
1645 errzDist=TMath::Sqrt(resol1*resol1+resol2*resol2);
1649 Double_t ex2 = pv->GetXRes();
1650 Double_t ey2= pv->GetYRes();
1651 errxDist=TMath::Sqrt(ex2*ex2+sigmax*sigmax);
1652 erryDist=TMath::Sqrt(ey2*ey2+sigmay*sigmay);
1654 if(nc2>=ncont && distanceZ>nSigmaDeltaZ*errzDist && distanceX<nSigmaXY*errxDist && distanceY<nSigmaXY*erryDist && distanceZ>distz)