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