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