1 /**************************************************************************
2 * Copyright(c) 1998-2007, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 //-------------------------------------------------------------------------
20 // Author: Markus Oldenburg, CERN
21 //-------------------------------------------------------------------------
26 #include <TFriendElement.h>
29 #include "AliAODEvent.h"
30 #include "AliAODHeader.h"
31 #include "AliAODTrack.h"
35 // definition of std AOD member names
36 const char* AliAODEvent::fAODListName[kAODListN] = {"header",
48 //______________________________________________________________________________
49 AliAODEvent::AliAODEvent() :
51 fAODObjects(new TList()),
65 // default constructor
68 //______________________________________________________________________________
69 AliAODEvent::AliAODEvent(const AliAODEvent& aod):
71 fAODObjects(new TList()),
72 fAODFolder(new TFolder()),
73 fHeader(new AliAODHeader(*aod.fHeader)),
74 fTracks(new TClonesArray(*aod.fTracks)),
75 fVertices(new TClonesArray(*aod.fVertices)),
76 fV0s(new TClonesArray(*aod.fV0s)),
77 fTracklets(new AliAODTracklets(*aod.fTracklets)),
78 fJets(new TClonesArray(*aod.fJets)),
79 fEmcalCells(new AliAODCaloCells(*aod.fEmcalCells)),
80 fPhosCells(new AliAODCaloCells(*aod.fPhosCells)),
81 fCaloClusters(new TClonesArray(*aod.fCaloClusters)),
82 fFmdClusters(new TClonesArray(*aod.fFmdClusters)),
83 fPmdClusters(new TClonesArray(*aod.fPmdClusters))
90 AddObject(fTracklets);
92 AddObject(fEmcalCells);
93 AddObject(fPhosCells);
94 AddObject(fCaloClusters);
95 AddObject(fFmdClusters);
96 AddObject(fPmdClusters);
101 //______________________________________________________________________________
102 AliAODEvent & AliAODEvent::operator=(const AliAODEvent& aod) {
104 // Assignment operator
106 if(&aod == this) return *this;
107 AliVEvent::operator=(aod);
109 fAODObjects = new TList();
110 fAODFolder = new TFolder();
111 fHeader = new AliAODHeader(*aod.fHeader);
112 fTracks = new TClonesArray(*aod.fTracks);
113 fVertices = new TClonesArray(*aod.fVertices);
114 fV0s = new TClonesArray(*aod.fV0s);
115 fTracklets = new AliAODTracklets(*aod.fTracklets);
116 fJets = new TClonesArray(*aod.fJets);
117 fEmcalCells = new AliAODCaloCells(*aod.fEmcalCells);
118 fPhosCells = new AliAODCaloCells(*aod.fPhosCells);
119 fCaloClusters = new TClonesArray(*aod.fCaloClusters);
120 fFmdClusters = new TClonesArray(*aod.fFmdClusters);
121 fPmdClusters = new TClonesArray(*aod.fPmdClusters);
123 fAODObjects = new TList();
127 AddObject(fVertices);
129 AddObject(fTracklets);
131 AddObject(fEmcalCells);
132 AddObject(fPhosCells);
133 AddObject(fCaloClusters);
134 AddObject(fFmdClusters);
135 AddObject(fPmdClusters);
141 //______________________________________________________________________________
142 AliAODEvent::~AliAODEvent()
149 //______________________________________________________________________________
150 void AliAODEvent::AddObject(TObject* obj)
152 // Add an object to the list of objects.
153 // Please be aware that in order to increase performance you should
154 // refrain from using TObjArrays (if possible). Use TClonesArrays, instead.
156 fAODObjects->AddLast(obj);
159 //______________________________________________________________________________
160 void AliAODEvent::RemoveObject(TObject* obj)
162 // Removes an object from the list of objects.
164 fAODObjects->Remove(obj);
167 //______________________________________________________________________________
168 TObject *AliAODEvent::FindListObject(const char *objName)
170 // Return the pointer to the object with the given name.
172 return fAODObjects->FindObject(objName);
175 //______________________________________________________________________________
176 void AliAODEvent::CreateStdContent()
178 // create the standard AOD content and set pointers
180 // create standard objects and add them to the TList of objects
181 AddObject(new AliAODHeader());
182 AddObject(new TClonesArray("AliAODTrack", 0));
183 AddObject(new TClonesArray("AliAODVertex", 0));
184 AddObject(new TClonesArray("AliAODv0", 0));
185 AddObject(new AliAODTracklets());
186 AddObject(new TClonesArray("AliAODJet", 0));
187 AddObject(new AliAODCaloCells());
188 AddObject(new AliAODCaloCells());
189 AddObject(new TClonesArray("AliAODCaloCluster", 0));
190 AddObject(new TClonesArray("AliAODFmdCluster", 0));
191 AddObject(new TClonesArray("AliAODPmdCluster", 0));
195 // read back pointers
201 //______________________________________________________________________________
202 void AliAODEvent::SetStdNames()
204 // introduce the standard naming
206 if(fAODObjects->GetEntries()==kAODListN){
207 for(int i = 0;i < fAODObjects->GetEntries();i++){
208 TObject *fObj = fAODObjects->At(i);
209 if(fObj->InheritsFrom("TNamed")){
210 ((TNamed*)fObj)->SetName(fAODListName[i]);
212 else if(fObj->InheritsFrom("TClonesArray")){
213 ((TClonesArray*)fObj)->SetName(fAODListName[i]);
218 printf("%s:%d SetStdNames() Wrong number of Std Entries \n",(char*)__FILE__,__LINE__);
222 void AliAODEvent::CreateStdFolders()
224 // Create the standard folder structure
225 fAODFolder = gROOT->GetRootFolder()->AddFolder("AOD", "AOD");
226 if(fAODObjects->GetEntries()==kAODListN){
227 for(int i = 0;i < fAODObjects->GetEntries();i++){
228 TObject *fObj = fAODObjects->At(i);
229 if(fObj->InheritsFrom("TClonesArray")){
230 fAODFolder->AddFolder(fAODListName[i], fAODListName[i], (TCollection*) fObj);
232 fAODFolder->AddFolder(fAODListName[i], fAODListName[i], 0);
237 printf("%s:%d CreateStdFolders() Wrong number of Std Entries \n",(char*)__FILE__,__LINE__);
241 //______________________________________________________________________________
242 void AliAODEvent::GetStdContent()
244 // set pointers for standard content
246 fHeader = (AliAODHeader*)fAODObjects->FindObject("header");
247 fTracks = (TClonesArray*)fAODObjects->FindObject("tracks");
248 fVertices = (TClonesArray*)fAODObjects->FindObject("vertices");
249 fV0s = (TClonesArray*)fAODObjects->FindObject("v0s");
250 fTracklets = (AliAODTracklets*)fAODObjects->FindObject("tracklets");
251 fJets = (TClonesArray*)fAODObjects->FindObject("jets");
252 fEmcalCells = (AliAODCaloCells*)fAODObjects->FindObject("emcalCells");
253 fPhosCells = (AliAODCaloCells*)fAODObjects->FindObject("phosCells");
254 fCaloClusters = (TClonesArray*)fAODObjects->FindObject("caloClusters");
255 fFmdClusters = (TClonesArray*)fAODObjects->FindObject("fmdClusters");
256 fPmdClusters = (TClonesArray*)fAODObjects->FindObject("pmdClusters");
259 //______________________________________________________________________________
260 void AliAODEvent::ResetStd(Int_t trkArrSize,
268 // deletes content of standard arrays and resets size
270 if (trkArrSize > fTracks->GetSize())
271 fTracks->Expand(trkArrSize);
274 if (vtxArrSize > fVertices->GetSize())
275 fVertices->Expand(vtxArrSize);
278 if (v0ArrSize > fV0s->GetSize())
279 fV0s->Expand(v0ArrSize);
282 if (jetSize > fJets->GetSize())
283 fJets->Expand(jetSize);
285 fCaloClusters->Delete();
286 if (caloClusSize > fCaloClusters->GetSize())
287 fCaloClusters->Expand(caloClusSize);
289 fFmdClusters->Delete();
290 if (fmdClusSize > fFmdClusters->GetSize())
291 fFmdClusters->Expand(fmdClusSize);
293 fPmdClusters->Delete();
294 if (pmdClusSize > fPmdClusters->GetSize())
295 fPmdClusters->Expand(pmdClusSize);
297 // Reset the tracklets
298 fTracklets->DeleteContainer();
299 fPhosCells->DeleteContainer();
300 fEmcalCells->DeleteContainer();
304 void AliAODEvent::ClearStd()
306 // clears the standard arrays
310 fTracklets ->DeleteContainer();
312 fEmcalCells ->DeleteContainer();
313 fPhosCells ->DeleteContainer();
314 fCaloClusters ->Clear();
315 fFmdClusters ->Clear();
316 fPmdClusters ->Clear();
319 //______________________________________________________________________________
320 Int_t AliAODEvent::GetMuonTracks(TRefArray *muonTracks) const
322 // fills the provided TRefArray with all found muon tracks
326 AliAODTrack *track = 0;
327 for (Int_t iTrack = 0; iTrack < GetNTracks(); iTrack++) {
328 if ((track = GetTrack(iTrack))->IsMuonTrack()) {
329 muonTracks->Add(track);
333 return muonTracks->GetEntriesFast();
337 void AliAODEvent::ReadFromTree(TTree *tree)
339 // connects aod event to tree
342 Printf("%s %d AliAODEvent::ReadFromTree() Zero Pointer to Tree \n",(char*)__FILE__,__LINE__);
346 if(!tree->GetTree())tree->LoadTree(0);
348 // Try to find AliAODEvent
349 AliAODEvent *aodEvent = 0;
350 aodEvent = (AliAODEvent*)tree->GetTree()->GetUserInfo()->FindObject("AliAODEvent");
352 // Check if already connected to tree
353 TList* connectedList = (TList*) (tree->GetUserInfo()->FindObject("AODObjectsConnectedToTree"));
355 // If connected use the connected list if objects
356 fAODObjects->Delete();
357 fAODObjects = connectedList;
362 // prevent a memory leak when reading back the TList
365 // create a new TList from the UserInfo TList...
366 // copy constructor does not work...
367 fAODObjects = (TList*)(aodEvent->GetList()->Clone());
368 fAODObjects->SetOwner(kFALSE);
369 if(fAODObjects->GetEntries()<kAODListN){
370 printf("%s %d AliAODEvent::ReadFromTree() TList contains less than the standard contents %d < %d \n",
371 (char*)__FILE__,__LINE__,fAODObjects->GetEntries(),kAODListN);
374 // Let's find out whether we have friends
375 TList* friendL = tree->GetListOfFriends();
380 while ((fe = (TFriendElement*)next())){
381 aodEvent = (AliAODEvent*)(fe->GetTree()->GetUserInfo()->FindObject("AliAODEvent"));
383 printf("No UserInfo on tree \n");
386 TList* objL = (TList*)(aodEvent->GetList()->Clone());
387 printf("Get list of object from tree %d !!\n", objL->GetEntries());
388 TIter nextobject(objL);
390 while(obj = nextobject())
392 printf("Adding object from friend %s !\n", obj->GetName());
393 fAODObjects->Add(obj);
394 } // object "branch" loop
400 // set the branch addresses
401 TIter next(fAODObjects);
403 while((el=(TNamed*)next())){
404 TString bname(el->GetName());
405 // check if branch exists under this Name
406 TBranch *br = tree->GetBranch(bname.Data());
408 tree->SetBranchAddress(bname.Data(),fAODObjects->GetObjectRef(el));
411 br = tree->GetBranch(Form("%s.",bname.Data()));
413 tree->SetBranchAddress(Form("%s.",bname.Data()),fAODObjects->GetObjectRef(el));
416 printf("%s %d AliAODEvent::ReadFromTree() No Branch found with Name %s. \n",
417 (char*)__FILE__,__LINE__,bname.Data());
422 // when reading back we are not owner of the list
423 // must not delete it
424 fAODObjects->SetOwner(kFALSE);
425 fAODObjects->SetName("AODObjectsConnectedToTree");
426 // we are not owner of the list objects
427 // must not delete it
428 tree->GetUserInfo()->Add(fAODObjects);
431 // we can't get the list from the user data, create standard content
432 // and set it by hand
434 TIter next(fAODObjects);
436 while((el=(TNamed*)next())){
437 TString bname(el->GetName());
438 tree->SetBranchAddress(bname.Data(),fAODObjects->GetObjectRef(el));
441 // when reading back we are not owner of the list
442 // must not delete it
443 fAODObjects->SetOwner(kFALSE);
447 //______________________________________________________________________________
448 void AliAODEvent::Print(Option_t *) const
450 // Something meaningful should be implemented here.