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