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