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