]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliESDEvent.cxx
New classes for event plane calculation (Johanna)
[u/mrichter/AliRoot.git] / STEER / AliESDEvent.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16 /* $Id$ */
17
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.
22 /*
23    AliESDEvent *ev= new AliESDEvent();
24    ev->ReadFromTree(esdTree);
25    ...
26     for (Int_t i=0; i<nev; i++) {
27       esdTree->GetEntry(i);
28       if(ev->GetAliESDOld())ev->CopyFromOldESD();
29 */
30 //   The AliESDInputHandler does this automatically for you
31 //
32 // Origin: Christian Klein-Boesing, CERN, Christian.Klein-Boesing@cern.ch
33 //-----------------------------------------------------------------
34
35 #include "TList.h"
36 #include "TRefArray.h"
37 #include <TNamed.h>
38 #include <TROOT.h>
39 #include <TInterpreter.h>
40
41 #include "AliESDEvent.h"
42 #include "AliESDfriend.h"
43 #include "AliESDVZERO.h"
44 #include "AliESDFMD.h"
45 #include "AliESD.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"
61 #include "AliESDv0.h"
62 #include "AliESDFMD.h"
63 #include "AliESDVZERO.h"
64 #include "AliMultiplicity.h"
65 #include "AliRawDataErrorLog.h"
66 #include "AliLog.h"
67 #include "AliESDACORDE.h"
68 #include "AliESDHLTDecision.h"
69 #include "AliCentrality.h"
70 #include "AliEventplane.h"
71
72 ClassImp(AliESDEvent)
73
74
75
76 // here we define the names, some classes are no TNamed, therefore the classnames 
77 // are the Names
78   const char* AliESDEvent::fgkESDListName[kESDListN] = {"AliESDRun",
79                                                         "AliESDHeader",
80                                                         "AliESDZDC",
81                                                         "AliESDFMD",
82                                                         "AliESDVZERO",
83                                                         "AliESDTZERO",
84                                                         "TPCVertex",
85                                                         "SPDVertex",
86                                                         "PrimaryVertex",
87                                                         "AliMultiplicity",
88                                                         "PHOSTrigger",
89                                                         "EMCALTrigger",
90                                                         "SPDPileupVertices",
91                                                         "TrkPileupVertices",
92                                                         "Tracks",
93                                                         "MuonTracks",
94                                                         "PmdTracks",
95                                                         "TrdTracks",
96                                                         "V0s",
97                                                         "Cascades",
98                                                         "Kinks",
99                                                         "CaloClusters",
100                                                         "EMCALCells",
101                                                         "PHOSCells",
102                                                         "AliRawDataErrorLogs",
103                                                         "AliESDACORDE",
104                                                         "AliTOFHeader"};
105
106 //______________________________________________________________________________
107 AliESDEvent::AliESDEvent():
108   AliVEvent(),
109   fESDObjects(new TList()),
110   fESDRun(0),
111   fHeader(0),
112   fESDZDC(0),
113   fESDFMD(0),
114   fESDVZERO(0),
115   fESDTZERO(0),
116   fTPCVertex(0),
117   fSPDVertex(0),
118   fPrimaryVertex(0),
119   fSPDMult(0),
120   fPHOSTrigger(0),
121   fEMCALTrigger(0),
122   fESDACORDE(0),
123   fSPDPileupVertices(0),
124   fTrkPileupVertices(0),
125   fTracks(0),
126   fMuonTracks(0),
127   fPmdTracks(0),
128   fTrdTracks(0),
129   fV0s(0),  
130   fCascades(0),
131   fKinks(0),
132   fCaloClusters(0),
133   fEMCALCells(0), fPHOSCells(0),
134   fErrorLogs(0),
135   fESDOld(0),
136   fESDFriendOld(0),
137   fConnected(kFALSE),
138   fUseOwnList(kFALSE),
139   fTOFHeader(0),
140   fCentrality(0),
141   fEventplane(0)
142 {
143 }
144 //______________________________________________________________________________
145 AliESDEvent::AliESDEvent(const AliESDEvent& esd):
146   AliVEvent(esd),
147   fESDObjects(new TList()),
148   fESDRun(new AliESDRun(*esd.fESDRun)),
149   fHeader(new AliESDHeader(*esd.fHeader)),
150   fESDZDC(new AliESDZDC(*esd.fESDZDC)),
151   fESDFMD(new AliESDFMD(*esd.fESDFMD)),
152   fESDVZERO(new AliESDVZERO(*esd.fESDVZERO)),
153   fESDTZERO(new AliESDTZERO(*esd.fESDTZERO)),
154   fTPCVertex(new AliESDVertex(*esd.fTPCVertex)),
155   fSPDVertex(new AliESDVertex(*esd.fSPDVertex)),
156   fPrimaryVertex(new AliESDVertex(*esd.fPrimaryVertex)),
157   fSPDMult(new AliMultiplicity(*esd.fSPDMult)),
158   fPHOSTrigger(new AliESDCaloTrigger(*esd.fPHOSTrigger)),
159   fEMCALTrigger(new AliESDCaloTrigger(*esd.fEMCALTrigger)),
160   fESDACORDE(new AliESDACORDE(*esd.fESDACORDE)),
161   fSPDPileupVertices(new TClonesArray(*esd.fSPDPileupVertices)),
162   fTrkPileupVertices(new TClonesArray(*esd.fTrkPileupVertices)),
163   fTracks(new TClonesArray(*esd.fTracks)),
164   fMuonTracks(new TClonesArray(*esd.fMuonTracks)),
165   fPmdTracks(new TClonesArray(*esd.fPmdTracks)),
166   fTrdTracks(new TClonesArray(*esd.fTrdTracks)),
167   fV0s(new TClonesArray(*esd.fV0s)),  
168   fCascades(new TClonesArray(*esd.fCascades)),
169   fKinks(new TClonesArray(*esd.fKinks)),
170   fCaloClusters(new TClonesArray(*esd.fCaloClusters)),
171   fEMCALCells(new AliESDCaloCells(*esd.fEMCALCells)),
172   fPHOSCells(new AliESDCaloCells(*esd.fPHOSCells)),
173   fErrorLogs(new TClonesArray(*esd.fErrorLogs)),
174   fESDOld(esd.fESDOld ? new AliESD(*esd.fESDOld) : 0),
175   fESDFriendOld(esd.fESDFriendOld ? new AliESDfriend(*esd.fESDFriendOld) : 0),
176   fConnected(esd.fConnected),
177   fUseOwnList(esd.fUseOwnList),
178   fTOFHeader(new AliTOFHeader(*esd.fTOFHeader)),
179   fCentrality(new AliCentrality(*esd.fCentrality)),
180   fEventplane(new AliEventplane(*esd.fEventplane))
181 {
182   // CKB init in the constructor list and only add here ...
183   AddObject(fESDRun);
184   AddObject(fHeader);
185   AddObject(fESDZDC);
186   AddObject(fESDFMD);
187   AddObject(fESDVZERO);
188   AddObject(fESDTZERO);
189   AddObject(fTPCVertex);
190   AddObject(fSPDVertex);
191   AddObject(fPrimaryVertex);
192   AddObject(fSPDMult);
193   AddObject(fPHOSTrigger);
194   AddObject(fEMCALTrigger);
195   AddObject(fSPDPileupVertices);
196   AddObject(fTrkPileupVertices);
197   AddObject(fTracks);
198   AddObject(fMuonTracks);
199   AddObject(fPmdTracks);
200   AddObject(fTrdTracks);
201   AddObject(fV0s);
202   AddObject(fCascades);
203   AddObject(fKinks);
204   AddObject(fCaloClusters);
205   AddObject(fEMCALCells);
206   AddObject(fPHOSCells);
207   AddObject(fErrorLogs);
208   AddObject(fESDACORDE);
209   AddObject(fTOFHeader);
210
211   GetStdContent();
212
213 }
214
215 //______________________________________________________________________________
216 AliESDEvent & AliESDEvent::operator=(const AliESDEvent& source) {
217
218   // Assignment operator
219
220   if(&source == this) return *this;
221   AliVEvent::operator=(source);
222
223   // This assumes that the list is already created
224   // and that the virtual void Copy(Tobject&) function
225   // is correctly implemented in the derived class
226   // otherwise only TObject::Copy() will be used
227
228
229
230   if((fESDObjects->GetSize()==0)&&(source.fESDObjects->GetSize()>=kESDListN)){
231     // We cover the case that we do not yet have the 
232     // standard content but the source has it
233     CreateStdContent();
234   }
235
236   TIter next(source.GetList());
237   TObject *its = 0;
238   TString name;
239   while ((its = next())) {
240     name.Form("%s", its->GetName());
241     TObject *mine = fESDObjects->FindObject(name.Data());
242     if(!mine){
243       TClass* pClass=TClass::GetClass(its->ClassName());
244       if (!pClass) {
245         AliWarning(Form("Can not find class description for entry %s (%s)\n",
246                         its->ClassName(), name.Data()));
247         continue;
248       }
249
250       mine=(TObject*)pClass->New();
251       if(!mine){
252       // not in this: can be added to list
253         AliWarning(Form("%s:%d Could not find %s for copying \n",
254                         (char*)__FILE__,__LINE__,name.Data()));
255         continue;
256       }  
257       if(mine->InheritsFrom("TNamed")){
258         ((TNamed*)mine)->SetName(name);
259       }
260       else if(mine->InheritsFrom("TCollection")){
261         if(mine->InheritsFrom("TClonesArray")) {
262           TClonesArray* tcits = dynamic_cast<TClonesArray*>(its);
263           if (tcits)
264             dynamic_cast<TClonesArray*>(mine)->SetClass(tcits->GetClass());
265         }
266         dynamic_cast<TCollection*>(mine)->SetName(name);
267       }
268       AliDebug(1, Form("adding object %s of type %s", mine->GetName(), mine->ClassName()));
269       AddObject(mine);
270     }  
271    
272     if(!its->InheritsFrom("TCollection")){
273       // simple objects
274       its->Copy(*mine);
275     }
276     else if(its->InheritsFrom("TClonesArray")){
277       // Create or expand the tclonesarray pointers
278       // so we can directly copy to the object
279       TClonesArray *its_tca = (TClonesArray*)its;
280       TClonesArray *mine_tca = (TClonesArray*)mine;
281
282       // this leaves the capacity of the TClonesArray the same
283       // except for a factor of 2 increase when size > capacity
284       // does not release any memory occupied by the tca
285       mine_tca->ExpandCreate(its_tca->GetEntriesFast());
286       for(int i = 0;i < its_tca->GetEntriesFast();++i){
287         // copy 
288         TObject *mine_tca_obj = mine_tca->At(i);
289         TObject *its_tca_obj = its_tca->At(i);
290         // no need to delete first
291         // pointers within the class should be handled by Copy()...
292         // Can there be Empty slots?
293         its_tca_obj->Copy(*mine_tca_obj);
294       }
295     }
296     else{
297       AliWarning(Form("%s:%d cannot copy TCollection \n",
298                       (char*)__FILE__,__LINE__));
299     }
300   }
301
302   fCentrality = source.fCentrality;
303   fEventplane = source.fEventplane;
304
305   fConnected  = source.fConnected;
306   fUseOwnList = source.fUseOwnList;
307
308   return *this;
309 }
310
311
312 //______________________________________________________________________________
313 AliESDEvent::~AliESDEvent()
314 {
315   //
316   // Standard destructor
317   //
318
319   // everthing on the list gets deleted automatically
320
321   
322   if(fESDObjects&&!fConnected)
323     {
324       delete fESDObjects;
325       fESDObjects = 0;
326     }
327   if (fCentrality) delete fCentrality;
328   if (fEventplane) delete fEventplane;
329   
330 }
331
332 void AliESDEvent::Copy(TObject &obj) const {
333
334   // interface to TOBject::Copy
335   // Copies the content of this into obj!
336   // bascially obj = *this
337
338   if(this==&obj)return;
339   AliESDEvent *robj = dynamic_cast<AliESDEvent*>(&obj);
340   if(!robj)return; // not an AliESEvent
341   *robj = *this;
342   return;
343 }
344
345 //______________________________________________________________________________
346 void AliESDEvent::Reset()
347 {
348
349   // Handle the cases
350   // Std content + Non std content
351
352   // Reset the standard contents
353   ResetStdContent(); 
354
355   //  reset for the old data without AliESDEvent...
356   if(fESDOld)fESDOld->Reset();
357   if(fESDFriendOld){
358     fESDFriendOld->~AliESDfriend();
359     new (fESDFriendOld) AliESDfriend();
360   }
361   // 
362
363   if(fESDObjects->GetSize()>kESDListN){
364     // we have non std content
365     // this also covers esdfriends
366     for(int i = kESDListN;i < fESDObjects->GetSize();++i){
367       TObject *pObject = fESDObjects->At(i);
368       // TClonesArrays
369       if(pObject->InheritsFrom(TClonesArray::Class())){
370         ((TClonesArray*)pObject)->Delete();
371       }
372       else if(!pObject->InheritsFrom(TCollection::Class())){
373         TClass *pClass = TClass::GetClass(pObject->ClassName());
374         if (pClass && pClass->GetListOfMethods()->FindObject("Clear")) {
375           AliDebug(1, Form("Clear for object %s class %s", pObject->GetName(), pObject->ClassName()));
376           pObject->Clear();
377         }
378         else {
379           AliDebug(1, Form("ResetWithPlacementNew for object %s class %s", pObject->GetName(), pObject->ClassName()));
380           ResetWithPlacementNew(pObject);
381         }
382       }
383       else{
384         AliWarning(Form("No reset for %s \n",
385                         pObject->ClassName()));
386       }
387     }
388   }
389
390 }
391
392 Bool_t AliESDEvent::ResetWithPlacementNew(TObject *pObject){
393   Long_t dtoronly = TObject::GetDtorOnly();
394   TClass *pClass = TClass::GetClass(pObject->ClassName()); 
395   TObject::SetDtorOnly(pObject);
396   delete pObject;
397   // Recreate with placement new
398   pClass->New(pObject);
399   // Restore the state.
400   TObject::SetDtorOnly((void*)dtoronly);
401   return kTRUE;
402 }
403
404 void AliESDEvent::ResetStdContent()
405 {
406   // Reset the standard contents
407   if(fESDRun) fESDRun->Reset();
408   if(fHeader) fHeader->Reset();
409   if(fESDZDC) fESDZDC->Reset();
410   if(fESDFMD) {
411     fESDFMD->Clear();
412   }
413   if(fESDVZERO){
414     // reset by callin d'to /c'tor keep the pointer
415     fESDVZERO->~AliESDVZERO();
416     new (fESDVZERO) AliESDVZERO();
417   }  
418   if(fESDACORDE){
419     fESDACORDE->~AliESDACORDE();
420     new (fESDACORDE) AliESDACORDE();    
421   } 
422   if(fESDTZERO) fESDTZERO->Reset(); 
423   // CKB no clear/reset implemented
424   if(fTPCVertex){
425     fTPCVertex->~AliESDVertex();
426     new (fTPCVertex) AliESDVertex();
427     fTPCVertex->SetName(fgkESDListName[kTPCVertex]);
428   }
429   if(fSPDVertex){
430     fSPDVertex->~AliESDVertex();
431     new (fSPDVertex) AliESDVertex();
432     fSPDVertex->SetName(fgkESDListName[kSPDVertex]);
433   }
434   if(fPrimaryVertex){
435     fPrimaryVertex->~AliESDVertex();
436     new (fPrimaryVertex) AliESDVertex();
437     fPrimaryVertex->SetName(fgkESDListName[kPrimaryVertex]);
438   }
439   if(fSPDMult){
440     fSPDMult->~AliMultiplicity();
441     new (fSPDMult) AliMultiplicity();
442   }
443   if(fTOFHeader){
444     fTOFHeader->~AliTOFHeader();
445     new (fTOFHeader) AliTOFHeader();
446     //fTOFHeader->SetName(fgkESDListName[kTOFHeader]);
447   }
448   if(fPHOSTrigger)fPHOSTrigger->DeAllocate(); 
449   if(fEMCALTrigger)fEMCALTrigger->DeAllocate(); 
450   if(fSPDPileupVertices)fSPDPileupVertices->Delete();
451   if(fTrkPileupVertices)fTrkPileupVertices->Delete();
452   if(fTracks)fTracks->Delete();
453   if(fMuonTracks)fMuonTracks->Delete();
454   if(fPmdTracks)fPmdTracks->Delete();
455   if(fTrdTracks)fTrdTracks->Delete();
456   if(fV0s)fV0s->Delete();
457   if(fCascades)fCascades->Delete();
458   if(fKinks)fKinks->Delete();
459   if(fCaloClusters)fCaloClusters->Delete();
460   if(fPHOSCells)fPHOSCells->DeleteContainer();
461   if(fEMCALCells)fEMCALCells->DeleteContainer();
462   if(fErrorLogs) fErrorLogs->Delete();
463
464   // don't reset fconnected fConnected and the list
465
466 }
467
468
469 Int_t AliESDEvent::AddV0(const AliESDv0 *v) {
470   //
471   // Add V0
472   //
473   TClonesArray &fv = *fV0s;
474   Int_t idx=fV0s->GetEntriesFast();
475   new(fv[idx]) AliESDv0(*v);
476   return idx;
477 }  
478
479 //______________________________________________________________________________
480 void AliESDEvent::Print(Option_t *) const 
481 {
482   //
483   // Print header information of the event
484   //
485   printf("ESD run information\n");
486   printf("Event # in file %d Bunch crossing # %d Orbit # %d Period # %d Run # %d Trigger %lld Magnetic field %f \n",
487          GetEventNumberInFile(),
488          GetBunchCrossNumber(),
489          GetOrbitNumber(),
490          GetPeriodNumber(),
491          GetRunNumber(),
492          GetTriggerMask(),
493          GetMagneticField() );
494   if (fPrimaryVertex)
495     printf("Vertex: (%.4f +- %.4f, %.4f +- %.4f, %.4f +- %.4f) cm\n",
496            fPrimaryVertex->GetXv(), fPrimaryVertex->GetXRes(),
497            fPrimaryVertex->GetYv(), fPrimaryVertex->GetYRes(),
498            fPrimaryVertex->GetZv(), fPrimaryVertex->GetZRes());
499   printf("Mean vertex in RUN: X=%.4f Y=%.4f Z=%.4f cm\n",
500          GetDiamondX(),GetDiamondY(),GetDiamondZ());
501   if(fSPDMult)
502     printf("SPD Multiplicity. Number of tracklets %d \n",
503            fSPDMult->GetNumberOfTracklets());
504   printf("Number of pileup primary vertices reconstructed with SPD %d\n", 
505          GetNumberOfPileupVerticesSPD());
506   printf("Number of pileup primary vertices reconstructed using the tracks %d\n",
507          GetNumberOfPileupVerticesTracks());
508   printf("Number of tracks: \n");
509   printf("                 charged   %d\n", GetNumberOfTracks());
510   printf("                 muon      %d\n", GetNumberOfMuonTracks());
511   printf("                 pmd       %d\n", GetNumberOfPmdTracks());
512   printf("                 trd       %d\n", GetNumberOfTrdTracks());
513   printf("                 v0        %d\n", GetNumberOfV0s());
514   printf("                 cascades  %d\n", GetNumberOfCascades());
515   printf("                 kinks     %d\n", GetNumberOfKinks());
516   if(fPHOSCells)printf("                 PHOSCells %d\n", fPHOSCells->GetNumberOfCells());
517   else printf("                 PHOSCells not in the Event\n");
518   if(fEMCALCells)printf("                 EMCALCells %d\n", fEMCALCells->GetNumberOfCells());
519   else printf("                 EMCALCells not in the Event\n");
520   printf("                 CaloClusters %d\n", GetNumberOfCaloClusters());
521   printf("                 FMD       %s\n", (fESDFMD ? "yes" : "no"));
522   printf("                 VZERO     %s\n", (fESDVZERO ? "yes" : "no"));
523   TObject* pHLTDecision=GetHLTTriggerDecision();
524   printf("HLT trigger decision: %s\n", pHLTDecision?pHLTDecision->GetOption():"not available");
525   if (pHLTDecision) pHLTDecision->Print("compact");
526
527   return;
528 }
529
530 void AliESDEvent::SetESDfriend(const AliESDfriend *ev) const {
531   //
532   // Attaches the complementary info to the ESD
533   //
534   if (!ev) return;
535
536   // to be sure that we set the tracks also
537   // in case of old esds 
538   // if(fESDOld)CopyFromOldESD();
539
540   Int_t ntrk=ev->GetNumberOfTracks();
541  
542   for (Int_t i=0; i<ntrk; i++) {
543     const AliESDfriendTrack *f=ev->GetTrack(i);
544     GetTrack(i)->SetFriendTrack(f);
545   }
546 }
547
548 Bool_t  AliESDEvent::RemoveKink(Int_t rm) const {
549   // ---------------------------------------------------------
550   // Remove a kink candidate and references to it from ESD,
551   // if this candidate does not come from a reconstructed decay
552   // Not yet implemented...
553   // ---------------------------------------------------------
554   Int_t last=GetNumberOfKinks()-1;
555   if ((rm<0)||(rm>last)) return kFALSE;
556
557   return kTRUE;
558 }
559
560 Bool_t  AliESDEvent::RemoveV0(Int_t rm) const {
561   // ---------------------------------------------------------
562   // Remove a V0 candidate and references to it from ESD,
563   // if this candidate does not come from a reconstructed decay
564   // ---------------------------------------------------------
565   Int_t last=GetNumberOfV0s()-1;
566   if ((rm<0)||(rm>last)) return kFALSE;
567
568   AliESDv0 *v0=GetV0(rm);
569   Int_t idxP=v0->GetPindex(), idxN=v0->GetNindex();
570
571   v0=GetV0(last);
572   Int_t lastIdxP=v0->GetPindex(), lastIdxN=v0->GetNindex();
573
574   Int_t used=0;
575
576   // Check if this V0 comes from a reconstructed decay
577   Int_t ncs=GetNumberOfCascades();
578   for (Int_t n=0; n<ncs; n++) {
579     AliESDcascade *cs=GetCascade(n);
580
581     Int_t csIdxP=cs->GetPindex();
582     Int_t csIdxN=cs->GetNindex();
583
584     if (idxP==csIdxP)
585        if (idxN==csIdxN) return kFALSE;
586
587     if (csIdxP==lastIdxP)
588        if (csIdxN==lastIdxN) used++;
589   }
590
591   //Replace the removed V0 with the last V0 
592   TClonesArray &a=*fV0s;
593   delete a.RemoveAt(rm);
594
595   if (rm==last) return kTRUE;
596
597   //v0 is pointing to the last V0 candidate... 
598   new (a[rm]) AliESDv0(*v0);
599   delete a.RemoveAt(last);
600
601   if (!used) return kTRUE;
602   
603
604   // Remap the indices of the daughters of reconstructed decays
605   for (Int_t n=0; n<ncs; n++) {
606     AliESDcascade *cs=GetCascade(n);
607
608
609     Int_t csIdxP=cs->GetPindex();
610     Int_t csIdxN=cs->GetNindex();
611
612     if (csIdxP==lastIdxP)
613       if (csIdxN==lastIdxN) {
614          cs->AliESDv0::SetIndex(1,idxP);
615          cs->AliESDv0::SetIndex(0,idxN);
616          used--;
617          if (!used) return kTRUE;
618       }
619   }
620
621   return kTRUE;
622 }
623
624 Bool_t  AliESDEvent::RemoveTrack(Int_t rm) const {
625   // ---------------------------------------------------------
626   // Remove a track and references to it from ESD,
627   // if this track does not come from a reconstructed decay
628   // ---------------------------------------------------------
629   Int_t last=GetNumberOfTracks()-1;
630   if ((rm<0)||(rm>last)) return kFALSE;
631
632   Int_t used=0;
633
634   // Check if this track comes from the reconstructed primary vertices
635   if (fTPCVertex && fTPCVertex->GetStatus()) {
636      UShort_t *primIdx=fTPCVertex->GetIndices();
637      Int_t n=fTPCVertex->GetNIndices();
638      while (n--) {
639        Int_t idx=Int_t(primIdx[n]);
640        if (rm==idx) return kFALSE;
641        if (idx==last) used++; 
642      }
643   }
644   if (fPrimaryVertex && fPrimaryVertex->GetStatus()) {
645      UShort_t *primIdx=fPrimaryVertex->GetIndices();
646      Int_t n=fPrimaryVertex->GetNIndices();
647      while (n--) {
648        Int_t idx=Int_t(primIdx[n]);
649        if (rm==idx) return kFALSE;
650        if (idx==last) used++; 
651      }
652   }
653   
654   // Check if this track comes from a reconstructed decay
655   Int_t nv0=GetNumberOfV0s();
656   for (Int_t n=0; n<nv0; n++) {
657     AliESDv0 *v0=GetV0(n);
658
659     Int_t idx=v0->GetNindex();
660     if (rm==idx) return kFALSE;
661     if (idx==last) used++;
662
663     idx=v0->GetPindex();
664     if (rm==idx) return kFALSE;
665     if (idx==last) used++;
666   }
667
668   Int_t ncs=GetNumberOfCascades();
669   for (Int_t n=0; n<ncs; n++) {
670     AliESDcascade *cs=GetCascade(n);
671
672     Int_t idx=cs->GetIndex();
673     if (rm==idx) return kFALSE;
674     if (idx==last) used++;
675
676     AliESDv0 *v0=cs;
677     idx=v0->GetNindex();
678     if (rm==idx) return kFALSE;
679     if (idx==last) used++;
680
681     idx=v0->GetPindex();
682     if (rm==idx) return kFALSE;
683     if (idx==last) used++;
684   }
685
686   Int_t nkn=GetNumberOfKinks();
687   for (Int_t n=0; n<nkn; n++) {
688     AliESDkink *kn=GetKink(n);
689
690     Int_t idx=kn->GetIndex(0);
691     if (rm==idx) return kFALSE;
692     if (idx==last) used++;
693
694     idx=kn->GetIndex(1);
695     if (rm==idx) return kFALSE;
696     if (idx==last) used++;
697   }
698
699   // Check if this track is associated with a CaloCluster
700   Int_t ncl=GetNumberOfCaloClusters();
701   for (Int_t n=0; n<ncl; n++) {
702     AliESDCaloCluster *cluster=GetCaloCluster(n);
703     TArrayI *arr=cluster->GetTracksMatched();
704     Int_t s=arr->GetSize();
705     while (s--) {
706       Int_t idx=arr->At(s);
707       if (rm==idx) return kFALSE;
708       if (idx==last) used++;     
709     }
710   }
711
712
713
714   //Replace the removed track with the last track 
715   TClonesArray &a=*fTracks;
716   delete a.RemoveAt(rm);
717
718   if (rm==last) return kTRUE;
719
720   AliESDtrack *t=GetTrack(last);
721   t->SetID(rm);
722   new (a[rm]) AliESDtrack(*t);
723   delete a.RemoveAt(last);
724
725
726   if (!used) return kTRUE;
727   
728
729   // Remap the indices of the tracks used for the primary vertex reconstruction
730   if (fTPCVertex && fTPCVertex->GetStatus()) {
731      UShort_t *primIdx=fTPCVertex->GetIndices();
732      Int_t n=fTPCVertex->GetNIndices();
733      while (n--) {
734        Int_t idx=Int_t(primIdx[n]);
735        if (idx==last) {
736           primIdx[n]=Short_t(rm); 
737           used--;
738           if (!used) return kTRUE;
739        }
740      }
741   }  
742   if (fPrimaryVertex && fPrimaryVertex->GetStatus()) {
743      UShort_t *primIdx=fPrimaryVertex->GetIndices();
744      Int_t n=fPrimaryVertex->GetNIndices();
745      while (n--) {
746        Int_t idx=Int_t(primIdx[n]);
747        if (idx==last) {
748           primIdx[n]=Short_t(rm); 
749           used--;
750           if (!used) return kTRUE;
751        }
752      }
753   }  
754
755   // Remap the indices of the daughters of reconstructed decays
756   for (Int_t n=0; n<nv0; n++) {
757     AliESDv0 *v0=GetV0(n);
758     if (v0->GetIndex(0)==last) {
759        v0->SetIndex(0,rm);
760        used--;
761        if (!used) return kTRUE;
762     }
763     if (v0->GetIndex(1)==last) {
764        v0->SetIndex(1,rm);
765        used--;
766        if (!used) return kTRUE;
767     }
768   }
769
770   for (Int_t n=0; n<ncs; n++) {
771     AliESDcascade *cs=GetCascade(n);
772     if (cs->GetIndex()==last) {
773        cs->SetIndex(rm);
774        used--;
775        if (!used) return kTRUE;
776     }
777     AliESDv0 *v0=cs;
778     if (v0->GetIndex(0)==last) {
779        v0->SetIndex(0,rm);
780        used--;
781        if (!used) return kTRUE;
782     }
783     if (v0->GetIndex(1)==last) {
784        v0->SetIndex(1,rm);
785        used--;
786        if (!used) return kTRUE;
787     }
788   }
789
790   for (Int_t n=0; n<nkn; n++) {
791     AliESDkink *kn=GetKink(n);
792     if (kn->GetIndex(0)==last) {
793        kn->SetIndex(rm,0);
794        used--;
795        if (!used) return kTRUE;
796     }
797     if (kn->GetIndex(1)==last) {
798        kn->SetIndex(rm,1);
799        used--;
800        if (!used) return kTRUE;
801     }
802   }
803
804   // Remap the indices of the tracks accosicated with CaloClusters
805   for (Int_t n=0; n<ncl; n++) {
806     AliESDCaloCluster *cluster=GetCaloCluster(n);
807     TArrayI *arr=cluster->GetTracksMatched();
808     Int_t s=arr->GetSize();
809     while (s--) {
810       Int_t idx=arr->At(s);
811       if (idx==last) {
812          arr->AddAt(rm,s);
813          used--; 
814          if (!used) return kTRUE;
815       }
816     }
817   }
818
819   return kTRUE;
820 }
821
822
823 Bool_t AliESDEvent::Clean(Float_t *cleanPars) {
824   //
825   // Remove the data which are not needed for the physics analysis.
826   //
827   // 1) Cleaning the V0 candidates
828   //    ---------------------------
829   //    If the cosine of the V0 pointing angle "csp" and 
830   //    the DCA between the daughter tracks "dca" does not satisfy 
831   //    the conditions 
832   //
833   //     csp > cleanPars[1] + dca/cleanPars[0]*(1.- cleanPars[1])
834   //
835   //    an attempt to remove this V0 candidate from ESD is made.
836   //
837   //    The V0 candidate gets removed if it does not belong to any 
838   //    recosntructed cascade decay
839   //
840   //    12.11.2007, optimal values: cleanPars[0]=0.5, cleanPars[1]=0.999
841   //
842   // 2) Cleaning the tracks
843   //    ----------------------
844   //    If track's transverse parameter is larger than cleanPars[2]
845   //                       OR
846   //    track's longitudinal parameter is larger than cleanPars[3]
847   //    an attempt to remove this track from ESD is made.
848   //
849   //    The track gets removed if it does not come 
850   //    from a reconstructed decay
851   //
852   Bool_t rc=kFALSE;
853
854   Float_t dcaMax=cleanPars[0];
855   Float_t cspMin=cleanPars[1];
856
857   Int_t nV0s=GetNumberOfV0s();
858   for (Int_t i=nV0s-1; i>=0; i--) {
859     AliESDv0 *v0=GetV0(i);
860
861     Float_t dca=v0->GetDcaV0Daughters();
862     Float_t csp=v0->GetV0CosineOfPointingAngle();
863     Float_t cspcut=cspMin + dca/dcaMax*(1.-cspMin);
864     if (csp > cspcut) continue;
865     if (RemoveV0(i)) rc=kTRUE;
866   }
867
868
869   Float_t dmax=cleanPars[2], zmax=cleanPars[3];
870
871   const AliESDVertex *vertex=GetPrimaryVertexSPD();
872   Bool_t vtxOK=vertex->GetStatus();
873   
874   Int_t nTracks=GetNumberOfTracks();
875   for (Int_t i=nTracks-1; i>=0; i--) {
876     AliESDtrack *track=GetTrack(i);
877     Float_t xy,z; track->GetImpactParameters(xy,z);
878     if ((TMath::Abs(xy) > dmax) || (vtxOK && (TMath::Abs(z) > zmax))) {
879       if (RemoveTrack(i)) rc=kTRUE;
880     }
881   }
882
883   return rc;
884 }
885
886 Char_t  AliESDEvent::AddPileupVertexSPD(const AliESDVertex *vtx) 
887 {
888     // Add a pileup primary vertex reconstructed with SPD
889     TClonesArray &ftr = *fSPDPileupVertices;
890     Char_t n=Char_t(ftr.GetEntriesFast());
891     AliESDVertex *vertex = new(ftr[n]) AliESDVertex(*vtx);
892     vertex->SetID(n);
893     return n;
894 }
895
896 Char_t  AliESDEvent::AddPileupVertexTracks(const AliESDVertex *vtx) 
897 {
898     // Add a pileup primary vertex reconstructed with SPD
899     TClonesArray &ftr = *fTrkPileupVertices;
900     Char_t n=Char_t(ftr.GetEntriesFast());
901     AliESDVertex *vertex = new(ftr[n]) AliESDVertex(*vtx);
902     vertex->SetID(n);
903     return n;
904 }
905
906 Int_t  AliESDEvent::AddTrack(const AliESDtrack *t) 
907 {
908     // Add track
909     TClonesArray &ftr = *fTracks;
910     AliESDtrack * track = new(ftr[fTracks->GetEntriesFast()])AliESDtrack(*t);
911     track->SetID(fTracks->GetEntriesFast()-1);
912     return  track->GetID();    
913 }
914
915  void AliESDEvent::AddMuonTrack(const AliESDMuonTrack *t) 
916 {
917     TClonesArray &fmu = *fMuonTracks;
918     new(fmu[fMuonTracks->GetEntriesFast()]) AliESDMuonTrack(*t);
919 }
920
921 void AliESDEvent::AddPmdTrack(const AliESDPmdTrack *t) 
922 {
923   TClonesArray &fpmd = *fPmdTracks;
924   new(fpmd[fPmdTracks->GetEntriesFast()]) AliESDPmdTrack(*t);
925 }
926
927 void AliESDEvent::AddTrdTrack(const AliESDTrdTrack *t) 
928 {
929   TClonesArray &ftrd = *fTrdTracks;
930   new(ftrd[fTrdTracks->GetEntriesFast()]) AliESDTrdTrack(*t);
931 }
932
933
934
935
936 Int_t AliESDEvent::AddKink(const AliESDkink *c) 
937 {
938     // Add kink
939     TClonesArray &fk = *fKinks;
940     AliESDkink * kink = new(fk[fKinks->GetEntriesFast()]) AliESDkink(*c);
941     kink->SetID(fKinks->GetEntriesFast()); // CKB different from the other imps..
942     return fKinks->GetEntriesFast()-1;
943 }
944
945
946 void AliESDEvent::AddCascade(const AliESDcascade *c) 
947 {
948   TClonesArray &fc = *fCascades;
949   new(fc[fCascades->GetEntriesFast()]) AliESDcascade(*c);
950 }
951
952
953 Int_t AliESDEvent::AddCaloCluster(const AliESDCaloCluster *c) 
954 {
955     // Add calocluster
956     TClonesArray &fc = *fCaloClusters;
957     AliESDCaloCluster *clus = new(fc[fCaloClusters->GetEntriesFast()]) AliESDCaloCluster(*c);
958     clus->SetID(fCaloClusters->GetEntriesFast()-1);
959     return fCaloClusters->GetEntriesFast()-1;
960   }
961
962
963 void  AliESDEvent::AddRawDataErrorLog(const AliRawDataErrorLog *log) const {
964   TClonesArray &errlogs = *fErrorLogs;
965   new(errlogs[errlogs.GetEntriesFast()])  AliRawDataErrorLog(*log);
966 }
967
968 void AliESDEvent::SetZDCData(AliESDZDC * obj)
969
970   // use already allocated space
971   if(fESDZDC)
972     *fESDZDC = *obj;
973 }
974
975 void  AliESDEvent::SetPrimaryVertexTPC(const AliESDVertex *vertex) 
976 {
977   // Set the TPC vertex
978   // use already allocated space
979   if(fTPCVertex){
980     *fTPCVertex = *vertex;
981     fTPCVertex->SetName(fgkESDListName[kTPCVertex]);
982   }
983 }
984
985 void  AliESDEvent::SetPrimaryVertexSPD(const AliESDVertex *vertex) 
986 {
987   // Set the SPD vertex
988   // use already allocated space
989   if(fSPDVertex){
990     *fSPDVertex = *vertex;
991     fSPDVertex->SetName(fgkESDListName[kSPDVertex]);
992   }
993 }
994
995 void  AliESDEvent::SetPrimaryVertexTracks(const AliESDVertex *vertex) 
996 {
997   // Set the primary vertex reconstructed using he ESD tracks.
998   // use already allocated space
999   if(fPrimaryVertex){
1000     *fPrimaryVertex = *vertex;
1001     fPrimaryVertex->SetName(fgkESDListName[kPrimaryVertex]);
1002   }
1003 }
1004
1005 const AliESDVertex * AliESDEvent::GetPrimaryVertex() const 
1006 {
1007   //
1008   // Get the "best" available reconstructed primary vertex.
1009   //
1010   if(fPrimaryVertex){
1011     if (fPrimaryVertex->GetStatus()) return fPrimaryVertex;
1012   }
1013   if(fSPDVertex){
1014     if (fSPDVertex->GetStatus()) return fSPDVertex;
1015   }
1016   if(fTPCVertex) return fTPCVertex;
1017   
1018   AliWarning("No primary vertex available. Returning the \"default\"...");
1019   return fSPDVertex;
1020 }
1021
1022 AliESDVertex * AliESDEvent::PrimaryVertexTracksUnconstrained() const 
1023 {
1024   //
1025   // Removes diamond constraint from fPrimaryVertex (reconstructed with tracks)
1026   // Returns a AliESDVertex which has to be deleted by the user
1027   //
1028   if(!fPrimaryVertex) {
1029     AliWarning("No primary vertex from tracks available.");
1030     return 0;
1031   }
1032   if(!fPrimaryVertex->GetStatus()) {
1033     AliWarning("No primary vertex from tracks available.");
1034     return 0;
1035   }
1036
1037   AliVertexerTracks vertexer(GetMagneticField());
1038   Float_t diamondxyz[3]={(Float_t)GetDiamondX(),(Float_t)GetDiamondY(),0.};
1039   Float_t diamondcovxy[3]; GetDiamondCovXY(diamondcovxy);
1040   Float_t diamondcov[6]={diamondcovxy[0],diamondcovxy[1],diamondcovxy[2],0.,0.,7.};
1041   AliESDVertex *vertex = 
1042     (AliESDVertex*)vertexer.RemoveConstraintFromVertex(fPrimaryVertex,diamondxyz,diamondcov);
1043
1044   return vertex;
1045 }
1046
1047 void AliESDEvent::SetMultiplicity(const AliMultiplicity *mul) 
1048 {
1049   // Set the SPD Multiplicity
1050   if(fSPDMult){
1051     *fSPDMult = *mul;
1052   }
1053 }
1054
1055
1056 void AliESDEvent::SetFMDData(AliESDFMD * obj) 
1057
1058   // use already allocated space
1059   if(fESDFMD){
1060     *fESDFMD = *obj;
1061   }
1062 }
1063
1064 void AliESDEvent::SetVZEROData(AliESDVZERO * obj)
1065
1066   // use already allocated space
1067   if(fESDVZERO)
1068     *fESDVZERO = *obj;
1069 }
1070
1071 void AliESDEvent::SetACORDEData(AliESDACORDE * obj)
1072 {
1073   if(fESDACORDE)
1074     *fESDACORDE = *obj;
1075 }
1076
1077
1078 void AliESDEvent::GetESDfriend(AliESDfriend *ev) const 
1079 {
1080   //
1081   // Extracts the complementary info from the ESD
1082   //
1083   if (!ev) return;
1084
1085   Int_t ntrk=GetNumberOfTracks();
1086
1087   for (Int_t i=0; i<ntrk; i++) {
1088     AliESDtrack *t=GetTrack(i);
1089     const AliESDfriendTrack *f=t->GetFriendTrack();
1090     ev->AddTrack(f);
1091
1092     t->ReleaseESDfriendTrack();// Not to have two copies of "friendTrack"
1093
1094   }
1095
1096   AliESDfriend *fr = (AliESDfriend*)(const_cast<AliESDEvent*>(this)->FindListObject("AliESDfriend"));
1097   if (fr) ev->SetVZEROfriend(fr->GetVZEROfriend());
1098 }
1099
1100 void AliESDEvent::AddObject(TObject* obj) 
1101 {
1102   // Add an object to the list of object.
1103   // Please be aware that in order to increase performance you should
1104   // refrain from using TObjArrays (if possible). Use TClonesArrays, instead.
1105   fESDObjects->SetOwner(kTRUE);
1106   fESDObjects->AddLast(obj);
1107 }
1108
1109
1110 void AliESDEvent::GetStdContent() 
1111 {
1112   // set pointers for standard content
1113   // get by name much safer and not a big overhead since not called very often
1114  
1115   fESDRun = (AliESDRun*)fESDObjects->FindObject(fgkESDListName[kESDRun]);
1116   fHeader = (AliESDHeader*)fESDObjects->FindObject(fgkESDListName[kHeader]);
1117   fESDZDC = (AliESDZDC*)fESDObjects->FindObject(fgkESDListName[kESDZDC]);
1118   fESDFMD = (AliESDFMD*)fESDObjects->FindObject(fgkESDListName[kESDFMD]);
1119   fESDVZERO = (AliESDVZERO*)fESDObjects->FindObject(fgkESDListName[kESDVZERO]);
1120   fESDTZERO = (AliESDTZERO*)fESDObjects->FindObject(fgkESDListName[kESDTZERO]);
1121   fTPCVertex = (AliESDVertex*)fESDObjects->FindObject(fgkESDListName[kTPCVertex]);
1122   fSPDVertex = (AliESDVertex*)fESDObjects->FindObject(fgkESDListName[kSPDVertex]);
1123   fPrimaryVertex = (AliESDVertex*)fESDObjects->FindObject(fgkESDListName[kPrimaryVertex]);
1124   fSPDMult =       (AliMultiplicity*)fESDObjects->FindObject(fgkESDListName[kSPDMult]);
1125   fPHOSTrigger = (AliESDCaloTrigger*)fESDObjects->FindObject(fgkESDListName[kPHOSTrigger]);
1126   fEMCALTrigger = (AliESDCaloTrigger*)fESDObjects->FindObject(fgkESDListName[kEMCALTrigger]);
1127   fSPDPileupVertices = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kSPDPileupVertices]);
1128   fTrkPileupVertices = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kTrkPileupVertices]);
1129   fTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kTracks]);
1130   fMuonTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kMuonTracks]);
1131   fPmdTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kPmdTracks]);
1132   fTrdTracks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kTrdTracks]);
1133   fV0s = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kV0s]);
1134   fCascades = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kCascades]);
1135   fKinks = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kKinks]);
1136   fCaloClusters = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kCaloClusters]);
1137   fEMCALCells = (AliESDCaloCells*)fESDObjects->FindObject(fgkESDListName[kEMCALCells]);
1138   fPHOSCells = (AliESDCaloCells*)fESDObjects->FindObject(fgkESDListName[kPHOSCells]);
1139   fErrorLogs = (TClonesArray*)fESDObjects->FindObject(fgkESDListName[kErrorLogs]);
1140   fESDACORDE = (AliESDACORDE*)fESDObjects->FindObject(fgkESDListName[kESDACORDE]);
1141   fTOFHeader = (AliTOFHeader*)fESDObjects->FindObject(fgkESDListName[kTOFHeader]);
1142 }
1143
1144 void AliESDEvent::SetStdNames(){
1145   // Set the names of the standard contents
1146   // 
1147   if(fESDObjects->GetEntries()>=kESDListN){
1148     for(int i = 0;i < fESDObjects->GetEntries() && i<kESDListN;i++){
1149       TObject *fObj = fESDObjects->At(i);
1150       if(fObj->InheritsFrom("TNamed")){
1151         ((TNamed*)fObj)->SetName(fgkESDListName[i]);
1152       }
1153       else if(fObj->InheritsFrom("TClonesArray")){
1154         ((TClonesArray*)fObj)->SetName(fgkESDListName[i]);
1155       }
1156     }
1157   }
1158   else{
1159      AliWarning("Std Entries missing");
1160   }
1161
1162
1163
1164 void AliESDEvent::CreateStdContent(Bool_t bUseThisList){
1165   fUseOwnList = bUseThisList;
1166   CreateStdContent();
1167 }
1168
1169 void AliESDEvent::CreateStdContent() 
1170 {
1171   // create the standard AOD content and set pointers
1172
1173   // create standard objects and add them to the TList of objects
1174   AddObject(new AliESDRun());
1175   AddObject(new AliESDHeader());
1176   AddObject(new AliESDZDC());
1177   AddObject(new AliESDFMD());
1178   AddObject(new AliESDVZERO());
1179   AddObject(new AliESDTZERO());
1180   AddObject(new AliESDVertex());
1181   AddObject(new AliESDVertex());
1182   AddObject(new AliESDVertex());
1183   AddObject(new AliMultiplicity());
1184   AddObject(new AliESDCaloTrigger());
1185   AddObject(new AliESDCaloTrigger());
1186   AddObject(new TClonesArray("AliESDVertex",0));
1187   AddObject(new TClonesArray("AliESDVertex",0));
1188   AddObject(new TClonesArray("AliESDtrack",0));
1189   AddObject(new TClonesArray("AliESDMuonTrack",0));
1190   AddObject(new TClonesArray("AliESDPmdTrack",0));
1191   AddObject(new TClonesArray("AliESDTrdTrack",0));
1192   AddObject(new TClonesArray("AliESDv0",0));
1193   AddObject(new TClonesArray("AliESDcascade",0));
1194   AddObject(new TClonesArray("AliESDkink",0));
1195   AddObject(new TClonesArray("AliESDCaloCluster",0));
1196   AddObject(new AliESDCaloCells());
1197   AddObject(new AliESDCaloCells());
1198   AddObject(new TClonesArray("AliRawDataErrorLog",0));
1199   AddObject(new AliESDACORDE()); 
1200   AddObject(new AliTOFHeader());
1201
1202   // check the order of the indices against enum...
1203
1204   // set names
1205   SetStdNames();
1206   // read back pointers
1207   GetStdContent();
1208 }
1209
1210 TObject* AliESDEvent::FindListObject(const char *name) const {
1211 //
1212 // Find object with name "name" in the list of branches
1213 //
1214   if(fESDObjects){
1215     return fESDObjects->FindObject(name);
1216   }
1217   return 0;
1218
1219
1220 Int_t AliESDEvent::GetPHOSClusters(TRefArray *clusters) const
1221 {
1222   // fills the provided TRefArray with all found phos clusters
1223   
1224   clusters->Clear();
1225   
1226   AliESDCaloCluster *cl = 0;
1227   for (Int_t i = 0; i < GetNumberOfCaloClusters(); i++) {
1228     
1229     if ( (cl = GetCaloCluster(i)) ) {
1230       if (cl->IsPHOS()){
1231         clusters->Add(cl);
1232         AliDebug(1,Form("IsPHOS cluster %d Size: %d \n",i,clusters->GetEntriesFast()));
1233       }
1234     }
1235   }
1236   return clusters->GetEntriesFast();
1237 }
1238
1239 Int_t AliESDEvent::GetEMCALClusters(TRefArray *clusters) const
1240 {
1241   // fills the provided TRefArray with all found emcal clusters
1242
1243   clusters->Clear();
1244
1245   AliESDCaloCluster *cl = 0;
1246   for (Int_t i = 0; i < GetNumberOfCaloClusters(); i++) {
1247
1248     if ( (cl = GetCaloCluster(i)) ) {
1249       if (cl->IsEMCAL()){
1250         clusters->Add(cl);
1251         AliDebug(1,Form("IsEMCAL cluster %d Size: %d \n",i,clusters->GetEntriesFast()));
1252       }
1253     }
1254   }
1255   return clusters->GetEntriesFast();
1256 }
1257
1258 void AliESDEvent::WriteToTree(TTree* tree) const {
1259   // Book the branches as in TTree::Branch(TCollection*)
1260   // but add a "." at the end of top level branches which are
1261   // not a TClonesArray
1262
1263
1264   TString branchname;
1265   TIter next(fESDObjects);
1266   const Int_t kSplitlevel = 99; // default value in TTree::Branch()
1267   const Int_t kBufsize = 32000; // default value in TTree::Branch()
1268   TObject *obj = 0;
1269
1270   while ((obj = next())) {
1271     branchname.Form("%s", obj->GetName());
1272     if(branchname.CompareTo("AliESDfriend")==0)branchname = "ESDfriend.";
1273     if ((kSplitlevel > 1) &&  !obj->InheritsFrom(TClonesArray::Class())) {
1274       if(!branchname.EndsWith("."))branchname += ".";
1275     }
1276     if (!tree->FindBranch(branchname)) {
1277       // For the custom streamer to be called splitlevel
1278       // has to be negative, only needed for HLT
1279       Int_t splitLevel = (TString(obj->ClassName()) == "AliHLTGlobalTriggerDecision") ? -1 : kSplitlevel - 1;
1280       tree->Bronch(branchname, obj->ClassName(), fESDObjects->GetObjectRef(obj),kBufsize, splitLevel);
1281     }
1282   }
1283 }
1284
1285
1286 void AliESDEvent::ReadFromTree(TTree *tree, Option_t* opt){
1287 //
1288 // Connect the ESDEvent to a tree
1289 //
1290   if(!tree){
1291     AliWarning("AliESDEvent::ReadFromTree() Zero Pointer to Tree \n");
1292     return;
1293   }
1294   // load the TTree
1295   if(!tree->GetTree())tree->LoadTree(0);
1296
1297   // if we find the "ESD" branch on the tree we do have the old structure
1298   if(tree->GetBranch("ESD")) {
1299     char ** address  = (char **)(tree->GetBranch("ESD")->GetAddress());
1300     // do we have the friend branch
1301     TBranch * esdFB = tree->GetBranch("ESDfriend.");
1302     char ** addressF = 0;
1303     if(esdFB)addressF = (char **)(esdFB->GetAddress());
1304     if (!address) {
1305       AliInfo("AliESDEvent::ReadFromTree() Reading old Tree");
1306       tree->SetBranchAddress("ESD",       &fESDOld);
1307       if(esdFB){
1308         tree->SetBranchAddress("ESDfriend.",&fESDFriendOld);
1309       }
1310     } else {
1311       AliInfo("AliESDEvent::ReadFromTree() Reading old Tree");
1312       AliInfo("Branch already connected. Using existing branch address.");
1313       fESDOld       = (AliESD*)       (*address);
1314       // addressF can still be 0, since branch needs to switched on
1315       if(addressF)fESDFriendOld = (AliESDfriend*) (*addressF);
1316     }
1317                                        
1318     //  have already connected the old ESD structure... ?
1319     // reuse also the pointer of the AlliESDEvent
1320     // otherwise create new ones
1321     TList* connectedList = (TList*) (tree->GetUserInfo()->FindObject("ESDObjectsConnectedToTree"));
1322   
1323     if(connectedList){
1324       // If connected use the connected list of objects
1325       if(fESDObjects!= connectedList){
1326         // protect when called twice 
1327         fESDObjects->Delete();
1328         fESDObjects = connectedList;
1329       }
1330       GetStdContent(); 
1331
1332       
1333       // The pointer to the friend changes when called twice via InitIO
1334       // since AliESDEvent is deleted
1335       TObject* oldf = FindListObject("AliESDfriend");
1336       TObject* newf = 0;
1337       if(addressF){
1338         newf = (TObject*)*addressF;
1339       }
1340       if(newf!=0&&oldf!=newf){
1341         // remove the old reference
1342         // Should we also delete it? Or is this handled in TTree I/O
1343         // since it is created by the first SetBranchAddress
1344         fESDObjects->Remove(oldf);
1345         // add the new one 
1346         fESDObjects->Add(newf);
1347       }
1348       
1349       fConnected = true;
1350       return;
1351     }
1352     // else...    
1353     CreateStdContent(); // create for copy
1354     // if we have the esdfriend add it, so we always can access it via the userinfo
1355     if(fESDFriendOld)AddObject(fESDFriendOld);
1356     // we are not owner of the list objects 
1357     // must not delete it
1358     fESDObjects->SetOwner(kTRUE);
1359     fESDObjects->SetName("ESDObjectsConnectedToTree");
1360     tree->GetUserInfo()->Add(fESDObjects);
1361     fConnected = true;
1362     return;
1363   }
1364   
1365
1366     delete fESDOld;
1367     fESDOld = 0;
1368   // Try to find AliESDEvent
1369   AliESDEvent *esdEvent = 0;
1370   esdEvent = (AliESDEvent*)tree->GetTree()->GetUserInfo()->FindObject("AliESDEvent");
1371   if(esdEvent){   
1372       // Check if already connected to tree
1373     esdEvent->Reset();
1374     TList* connectedList = (TList*) (tree->GetUserInfo()->FindObject("ESDObjectsConnectedToTree"));
1375
1376     
1377     if (connectedList && (strcmp(opt, "reconnect"))) {
1378       // If connected use the connected list if objects
1379       fESDObjects->Delete();
1380       fESDObjects = connectedList;
1381       GetStdContent(); 
1382       fConnected = true;
1383       return;
1384     }
1385
1386     // Connect to tree
1387     // prevent a memory leak when reading back the TList
1388     // if (!(strcmp(opt, "reconnect"))) fESDObjects->Delete();
1389     
1390     if(!fUseOwnList){
1391       // create a new TList from the UserInfo TList... 
1392       // copy constructor does not work...
1393       fESDObjects = (TList*)(esdEvent->GetList()->Clone());
1394       fESDObjects->SetOwner(kTRUE);
1395     }
1396     else if ( fESDObjects->GetEntries()==0){
1397       // at least create the std content if we want to read to our list
1398       CreateStdContent(); 
1399     }
1400
1401     // in principle
1402     // we only need new things in the list if we do no already have it..
1403     // TODO just add new entries
1404
1405     if(fESDObjects->GetEntries()<kESDListN){
1406       AliWarning(Form("AliESDEvent::ReadFromTree() TList contains less than the standard contents %d < %d \n",
1407                       fESDObjects->GetEntries(),kESDListN));
1408     }
1409     // set the branch addresses
1410     TIter next(fESDObjects);
1411     TNamed *el;
1412     while((el=(TNamed*)next())){
1413       TString bname(el->GetName());
1414       if(bname.CompareTo("AliESDfriend")==0)
1415         {
1416           // AliESDfriend does not have a name ...
1417             TBranch *br = tree->GetBranch("ESDfriend.");
1418             if (br) tree->SetBranchAddress("ESDfriend.",fESDObjects->GetObjectRef(el));
1419         }
1420       else{
1421         // check if branch exists under this Name
1422         TBranch *br = tree->GetBranch(bname.Data());
1423         if(br){
1424           tree->SetBranchAddress(bname.Data(),fESDObjects->GetObjectRef(el));
1425         }
1426         else{
1427           br = tree->GetBranch(Form("%s.",bname.Data()));
1428           if(br){
1429             tree->SetBranchAddress(Form("%s.",bname.Data()),fESDObjects->GetObjectRef(el));
1430           }
1431           else{
1432             AliWarning(Form("AliESDEvent::ReadFromTree() No Branch found with Name %s or %s.",bname.Data(),bname.Data()));
1433           }
1434
1435         }
1436       }
1437     }
1438     GetStdContent();
1439     // when reading back we are not owner of the list 
1440     // must not delete it
1441     fESDObjects->SetOwner(kTRUE);
1442     fESDObjects->SetName("ESDObjectsConnectedToTree");
1443     // we are not owner of the list objects 
1444     // must not delete it
1445     tree->GetUserInfo()->Add(fESDObjects);
1446     tree->GetUserInfo()->SetOwner(kFALSE);
1447     fConnected = true;
1448   }// no esdEvent -->
1449   else {
1450     // we can't get the list from the user data, create standard content
1451     // and set it by hand (no ESDfriend at the moment
1452     CreateStdContent();
1453     TIter next(fESDObjects);
1454     TNamed *el;
1455     while((el=(TNamed*)next())){
1456       TString bname(el->GetName());    
1457       TBranch *br = tree->GetBranch(bname.Data());
1458       if(br){
1459         tree->SetBranchAddress(bname.Data(),fESDObjects->GetObjectRef(el));
1460       }
1461       else{
1462         br = tree->GetBranch(Form("%s.",bname.Data()));
1463         if(br){
1464           tree->SetBranchAddress(Form("%s.",bname.Data()),fESDObjects->GetObjectRef(el));
1465         }
1466       }
1467     }
1468     GetStdContent();
1469     // when reading back we are not owner of the list 
1470     // must not delete it
1471     fESDObjects->SetOwner(kTRUE);
1472   }
1473 }
1474
1475
1476 void AliESDEvent::CopyFromOldESD()
1477 {
1478   // Method which copies over everthing from the old esd structure to the 
1479   // new  
1480   if(fESDOld){
1481     ResetStdContent();
1482      // Run
1483     SetRunNumber(fESDOld->GetRunNumber());
1484     SetPeriodNumber(fESDOld->GetPeriodNumber());
1485     SetMagneticField(fESDOld->GetMagneticField());
1486   
1487     // leave out diamond ...
1488     // SetDiamond(const AliESDVertex *vertex) { fESDRun->SetDiamond(vertex);}
1489
1490     // header
1491     SetTriggerMask(fESDOld->GetTriggerMask());
1492     SetOrbitNumber(fESDOld->GetOrbitNumber());
1493     SetTimeStamp(fESDOld->GetTimeStamp());
1494     SetEventType(fESDOld->GetEventType());
1495     SetEventNumberInFile(fESDOld->GetEventNumberInFile());
1496     SetBunchCrossNumber(fESDOld->GetBunchCrossNumber());
1497     SetTriggerCluster(fESDOld->GetTriggerCluster());
1498
1499     // ZDC
1500
1501     SetZDC(fESDOld->GetZDCN1Energy(),
1502            fESDOld->GetZDCP1Energy(),
1503            fESDOld->GetZDCEMEnergy(),
1504            0,
1505            fESDOld->GetZDCN2Energy(),
1506            fESDOld->GetZDCP2Energy(),
1507            fESDOld->GetZDCParticipants(),
1508            0,
1509            0,
1510            0,
1511            0,
1512            0,
1513            0);
1514
1515     // FMD
1516     
1517     if(fESDOld->GetFMDData())SetFMDData(fESDOld->GetFMDData());
1518
1519     // T0
1520
1521     SetT0zVertex(fESDOld->GetT0zVertex());
1522     SetT0(fESDOld->GetT0());
1523     //  leave amps out
1524
1525     // VZERO
1526     if (fESDOld->GetVZEROData()) SetVZEROData(fESDOld->GetVZEROData());
1527
1528     if(fESDOld->GetVertex())SetPrimaryVertexSPD(fESDOld->GetVertex());
1529
1530     if(fESDOld->GetPrimaryVertex())SetPrimaryVertexTracks(fESDOld->GetPrimaryVertex());
1531
1532     if(fESDOld->GetMultiplicity())SetMultiplicity(fESDOld->GetMultiplicity());
1533
1534     for(int i = 0;i<fESDOld->GetNumberOfTracks();i++){
1535       AddTrack(fESDOld->GetTrack(i));
1536     }
1537
1538     for(int i = 0;i<fESDOld->GetNumberOfMuonTracks();i++){
1539       AddMuonTrack(fESDOld->GetMuonTrack(i));
1540     }
1541
1542     for(int i = 0;i<fESDOld->GetNumberOfPmdTracks();i++){
1543       AddPmdTrack(fESDOld->GetPmdTrack(i));
1544     }
1545
1546     for(int i = 0;i<fESDOld->GetNumberOfTrdTracks();i++){
1547       AddTrdTrack(fESDOld->GetTrdTrack(i));
1548     }
1549
1550     for(int i = 0;i<fESDOld->GetNumberOfV0s();i++){
1551       AddV0(fESDOld->GetV0(i));
1552     }
1553
1554     for(int i = 0;i<fESDOld->GetNumberOfCascades();i++){
1555       AddCascade(fESDOld->GetCascade(i));
1556     }
1557
1558     for(int i = 0;i<fESDOld->GetNumberOfKinks();i++){
1559       AddKink(fESDOld->GetKink(i));
1560     }
1561
1562
1563     for(int i = 0;i<fESDOld->GetNumberOfCaloClusters();i++){
1564       AddCaloCluster(fESDOld->GetCaloCluster(i));
1565     }
1566
1567   }// if fesdold
1568 }
1569
1570 Bool_t AliESDEvent::IsEventSelected(const char *trigExpr) const
1571 {
1572   // Check if the event satisfies the trigger
1573   // selection expression trigExpr.
1574   // trigExpr can be any logical expression
1575   // of the trigger classes defined in AliESDRun
1576   // In case of wrong syntax return kTRUE.
1577
1578   TString expr(trigExpr);
1579   if (expr.IsNull()) return kTRUE;
1580
1581   ULong64_t mask = GetTriggerMask();
1582   for(Int_t itrig = 0; itrig < AliESDRun::kNTriggerClasses; itrig++) {
1583     if (mask & (1ull << itrig)) {
1584       expr.ReplaceAll(GetESDRun()->GetTriggerClass(itrig),"1");
1585     }
1586     else {
1587       expr.ReplaceAll(GetESDRun()->GetTriggerClass(itrig),"0");
1588     }
1589   }
1590
1591   Int_t error;
1592   if ((gROOT->ProcessLineFast(expr.Data(),&error) == 0) &&
1593       (error == TInterpreter::kNoError)) {
1594     return kFALSE;
1595   }
1596
1597   return kTRUE;
1598
1599 }
1600
1601 TObject*  AliESDEvent::GetHLTTriggerDecision() const
1602 {
1603   // get the HLT trigger decission object
1604
1605   // cast away const'nes because the FindListObject method
1606   // is not const
1607   AliESDEvent* pNonConst=const_cast<AliESDEvent*>(this);
1608   return pNonConst->FindListObject("HLTGlobalTrigger");
1609 }
1610
1611 TString   AliESDEvent::GetHLTTriggerDescription() const
1612 {
1613   // get the HLT trigger decission description
1614   TString description;
1615   TObject* pDecision=GetHLTTriggerDecision();
1616   if (pDecision) {
1617     description=pDecision->GetTitle();
1618   }
1619
1620   return description;
1621 }
1622
1623 Bool_t    AliESDEvent::IsHLTTriggerFired(const char* name) const
1624 {
1625   // get the HLT trigger decission description
1626   TObject* pDecision=GetHLTTriggerDecision();
1627   if (!pDecision) return kFALSE;
1628
1629   Option_t* option=pDecision->GetOption();
1630   if (option==NULL || *option!='1') return kFALSE;
1631
1632   if (name) {
1633     TString description=GetHLTTriggerDescription();
1634     Int_t index=description.Index(name);
1635     if (index<0) return kFALSE;
1636     index+=strlen(name);
1637     if (index>=description.Length()) return kFALSE;
1638     if (description[index]!=0 && description[index]!=' ') return kFALSE;
1639   }
1640   return kTRUE;
1641 }
1642
1643 //______________________________________________________________________________
1644 Bool_t  AliESDEvent::IsPileupFromSPD(Int_t minContributors, 
1645                                      Double_t minZdist, 
1646                                      Double_t nSigmaZdist, 
1647                                      Double_t nSigmaDiamXY, 
1648                                      Double_t nSigmaDiamZ) const{
1649   //
1650   // This function checks if there was a pile up
1651   // reconstructed with SPD
1652   //
1653   Int_t nc1=fSPDVertex->GetNContributors();
1654   if(nc1<1) return kFALSE;
1655   Int_t nPileVert=GetNumberOfPileupVerticesSPD();
1656   if(nPileVert==0) return kFALSE;
1657   
1658   for(Int_t i=0; i<nPileVert;i++){
1659     const AliESDVertex* pv=GetPileupVertexSPD(i);
1660     Int_t nc2=pv->GetNContributors();
1661     if(nc2>=minContributors){
1662       Double_t z1=fSPDVertex->GetZ();
1663       Double_t z2=pv->GetZ();
1664       Double_t distZ=TMath::Abs(z2-z1);
1665       Double_t distZdiam=TMath::Abs(z2-GetDiamondZ());
1666       Double_t cutZdiam=nSigmaDiamZ*TMath::Sqrt(GetSigma2DiamondZ());
1667       if(GetSigma2DiamondZ()<0.0001)cutZdiam=99999.; //protection for missing z diamond information
1668       if(distZ>minZdist && distZdiam<cutZdiam){
1669         Double_t x2=pv->GetX();
1670         Double_t y2=pv->GetY();
1671         Double_t distXdiam=TMath::Abs(x2-GetDiamondX());
1672         Double_t distYdiam=TMath::Abs(y2-GetDiamondY());
1673         Double_t cov1[6],cov2[6];       
1674         fSPDVertex->GetCovarianceMatrix(cov1);
1675         pv->GetCovarianceMatrix(cov2);
1676         Double_t errxDist=TMath::Sqrt(cov2[0]+GetSigma2DiamondX());
1677         Double_t erryDist=TMath::Sqrt(cov2[2]+GetSigma2DiamondY());
1678         Double_t errzDist=TMath::Sqrt(cov1[5]+cov2[5]);
1679         Double_t cutXdiam=nSigmaDiamXY*errxDist;
1680         if(GetSigma2DiamondX()<0.0001)cutXdiam=99999.; //protection for missing diamond information
1681         Double_t cutYdiam=nSigmaDiamXY*erryDist;
1682         if(GetSigma2DiamondY()<0.0001)cutYdiam=99999.; //protection for missing diamond information
1683         if( (distXdiam<cutXdiam) && (distYdiam<cutYdiam) && (distZ>nSigmaZdist*errzDist) ){
1684           return kTRUE;
1685         }
1686       }
1687     }
1688   }
1689   return kFALSE;
1690 }
1691
1692 //______________________________________________________________________________
1693 void AliESDEvent::EstimateMultiplicity(Int_t &tracklets, Int_t &trITSTPC, Int_t &trITSSApure, Double_t eta, Bool_t useDCAFlag,Bool_t useV0Flag) const
1694 {
1695   //
1696   // calculates 3 estimators for the multiplicity in the -eta:eta range
1697   // tracklets   : using SPD tracklets only
1698   // trITSTPC    : using TPC/ITS + complementary ITS SA tracks + tracklets from clusters not used by tracks
1699   // trITSSApure : using ITS standalone tracks + tracklets from clusters not used by tracks
1700   // if useDCAFlag is true: account for the ESDtrack flag marking the tracks with large DCA
1701   // if useV0Flag  is true: account for the ESDtrack flag marking conversion and K0's V0s
1702   tracklets = trITSSApure = trITSTPC = 0;
1703   int ntr = fSPDMult ? fSPDMult->GetNumberOfTracklets() : 0;
1704   //
1705   // count tracklets
1706   for (int itr=ntr;itr--;) { 
1707     if (TMath::Abs(fSPDMult->GetEta(itr))>eta) continue;
1708     tracklets++;
1709     if (fSPDMult->FreeClustersTracklet(itr,0)) trITSTPC++;    // not used in ITS/TPC or ITS_SA track
1710     if (fSPDMult->FreeClustersTracklet(itr,1)) trITSSApure++; // not used in ITS_SA_Pure track
1711   }
1712   //
1713   // count real tracks
1714   ntr = GetNumberOfTracks();
1715   for (int itr=ntr;itr--;) {
1716     AliESDtrack *t = GetTrack(itr);
1717     if (TMath::Abs(t->Eta())>eta) continue;
1718     if (!t->IsOn(AliESDtrack::kITSin)) continue;
1719     if (useDCAFlag && t->IsOn(AliESDtrack::kMultSec))  continue;
1720     if (useV0Flag  && t->IsOn(AliESDtrack::kMultInV0)) continue;    
1721     if (t->IsOn(AliESDtrack::kITSpureSA)) trITSSApure++;
1722     else                                  trITSTPC++;
1723   }
1724   //
1725 }
1726
1727 Bool_t AliESDEvent::IsPileupFromSPDInMultBins() const {
1728     Int_t nTracklets=GetMultiplicity()->GetNumberOfTracklets();
1729     if(nTracklets<20) return IsPileupFromSPD(3,0.8);
1730     else if(nTracklets<50) return IsPileupFromSPD(4,0.8);
1731     else return IsPileupFromSPD(5,0.8);
1732 }
1733
1734 void  AliESDEvent::SetTOFHeader(const AliTOFHeader *header)
1735 {
1736   //
1737   // Set the TOF event_time
1738   //
1739
1740   if (fTOFHeader) {
1741     *fTOFHeader=*header;
1742     //fTOFHeader->SetName(fgkESDListName[kTOFHeader]);
1743   }
1744   else {
1745     // for analysis of reconstructed events
1746     // when this information is not avaliable
1747     fTOFHeader = new AliTOFHeader(*header);
1748     //AddObject(fTOFHeader);
1749   }
1750
1751 }
1752
1753 AliCentrality* AliESDEvent::GetCentrality()
1754 {
1755     if (!fCentrality) fCentrality = new AliCentrality();
1756     return  fCentrality;
1757 }
1758
1759 AliEventplane* AliESDEvent::GetEventplane()
1760 {
1761     if (!fEventplane) fEventplane = new AliEventplane();
1762     return  fEventplane;
1763 }