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 //-------------------------------------------------------------------------
27 #include "AliAODEvent.h"
28 #include "AliAODHeader.h"
29 #include "AliAODTrack.h"
33 // definition of std AOD member names
34 const char* AliAODEvent::fAODListName[kAODListN] = {"header",
46 //______________________________________________________________________________
47 AliAODEvent::AliAODEvent() :
49 fAODObjects(new TList()),
63 // default constructor
66 //______________________________________________________________________________
67 AliAODEvent::AliAODEvent(const AliAODEvent& aod):
69 fAODObjects(new TList()),
70 fAODFolder(new TFolder()),
71 fHeader(new AliAODHeader(*aod.fHeader)),
72 fTracks(new TClonesArray(*aod.fTracks)),
73 fVertices(new TClonesArray(*aod.fVertices)),
74 fV0s(new TClonesArray(*aod.fV0s)),
75 fTracklets(new AliAODTracklets(*aod.fTracklets)),
76 fJets(new TClonesArray(*aod.fJets)),
77 fEmcalCells(new AliAODCaloCells(*aod.fEmcalCells)),
78 fPhosCells(new AliAODCaloCells(*aod.fPhosCells)),
79 fCaloClusters(new TClonesArray(*aod.fCaloClusters)),
80 fFmdClusters(new TClonesArray(*aod.fFmdClusters)),
81 fPmdClusters(new TClonesArray(*aod.fPmdClusters))
88 AddObject(fTracklets);
90 AddObject(fEmcalCells);
91 AddObject(fPhosCells);
92 AddObject(fCaloClusters);
93 AddObject(fFmdClusters);
94 AddObject(fPmdClusters);
99 //______________________________________________________________________________
100 AliAODEvent & AliAODEvent::operator=(const AliAODEvent& aod) {
102 // Assignment operator
104 if(&aod == this) return *this;
105 AliVEvent::operator=(aod);
107 fAODObjects = new TList();
108 fAODFolder = new TFolder();
109 fHeader = new AliAODHeader(*aod.fHeader);
110 fTracks = new TClonesArray(*aod.fTracks);
111 fVertices = new TClonesArray(*aod.fVertices);
112 fV0s = new TClonesArray(*aod.fV0s);
113 fTracklets = new AliAODTracklets(*aod.fTracklets);
114 fJets = new TClonesArray(*aod.fJets);
115 fEmcalCells = new AliAODCaloCells(*aod.fEmcalCells);
116 fPhosCells = new AliAODCaloCells(*aod.fPhosCells);
117 fCaloClusters = new TClonesArray(*aod.fCaloClusters);
118 fFmdClusters = new TClonesArray(*aod.fFmdClusters);
119 fPmdClusters = new TClonesArray(*aod.fPmdClusters);
121 fAODObjects = new TList();
125 AddObject(fVertices);
127 AddObject(fTracklets);
129 AddObject(fEmcalCells);
130 AddObject(fPhosCells);
131 AddObject(fCaloClusters);
132 AddObject(fFmdClusters);
133 AddObject(fPmdClusters);
139 //______________________________________________________________________________
140 AliAODEvent::~AliAODEvent()
147 //______________________________________________________________________________
148 void AliAODEvent::AddObject(TObject* obj)
150 // Add an object to the list of objects.
151 // Please be aware that in order to increase performance you should
152 // refrain from using TObjArrays (if possible). Use TClonesArrays, instead.
154 fAODObjects->AddLast(obj);
157 //______________________________________________________________________________
158 void AliAODEvent::RemoveObject(TObject* obj)
160 // Removes an object from the list of objects.
162 fAODObjects->Remove(obj);
165 //______________________________________________________________________________
166 TObject *AliAODEvent::FindListObject(const char *objName)
168 // Return the pointer to the object with the given name.
170 return fAODObjects->FindObject(objName);
173 //______________________________________________________________________________
174 void AliAODEvent::CreateStdContent()
176 // create the standard AOD content and set pointers
178 // create standard objects and add them to the TList of objects
179 AddObject(new AliAODHeader());
180 AddObject(new TClonesArray("AliAODTrack", 0));
181 AddObject(new TClonesArray("AliAODVertex", 0));
182 AddObject(new TClonesArray("AliAODv0", 0));
183 AddObject(new AliAODTracklets());
184 AddObject(new TClonesArray("AliAODJet", 0));
185 AddObject(new AliAODCaloCells());
186 AddObject(new AliAODCaloCells());
187 AddObject(new TClonesArray("AliAODCaloCluster", 0));
188 AddObject(new TClonesArray("AliAODFmdCluster", 0));
189 AddObject(new TClonesArray("AliAODPmdCluster", 0));
193 // read back pointers
199 //______________________________________________________________________________
200 void AliAODEvent::SetStdNames()
202 // introduce the standard naming
204 if(fAODObjects->GetEntries()==kAODListN){
205 for(int i = 0;i < fAODObjects->GetEntries();i++){
206 TObject *fObj = fAODObjects->At(i);
207 if(fObj->InheritsFrom("TNamed")){
208 ((TNamed*)fObj)->SetName(fAODListName[i]);
210 else if(fObj->InheritsFrom("TClonesArray")){
211 ((TClonesArray*)fObj)->SetName(fAODListName[i]);
216 printf("%s:%d SetStdNames() Wrong number of Std Entries \n",(char*)__FILE__,__LINE__);
220 void AliAODEvent::CreateStdFolders()
222 // Create the standard folder structure
223 fAODFolder = gROOT->GetRootFolder()->AddFolder("AOD", "AOD");
224 if(fAODObjects->GetEntries()==kAODListN){
225 for(int i = 0;i < fAODObjects->GetEntries();i++){
226 TObject *fObj = fAODObjects->At(i);
227 if(fObj->InheritsFrom("TClonesArray")){
228 fAODFolder->AddFolder(fAODListName[i], fAODListName[i], (TCollection*) fObj);
230 fAODFolder->AddFolder(fAODListName[i], fAODListName[i], 0);
235 printf("%s:%d CreateStdFolders() Wrong number of Std Entries \n",(char*)__FILE__,__LINE__);
239 //______________________________________________________________________________
240 void AliAODEvent::GetStdContent()
242 // set pointers for standard content
244 fHeader = (AliAODHeader*)fAODObjects->FindObject("header");
245 fTracks = (TClonesArray*)fAODObjects->FindObject("tracks");
246 fVertices = (TClonesArray*)fAODObjects->FindObject("vertices");
247 fV0s = (TClonesArray*)fAODObjects->FindObject("v0s");
248 fTracklets = (AliAODTracklets*)fAODObjects->FindObject("tracklets");
249 fJets = (TClonesArray*)fAODObjects->FindObject("jets");
250 fEmcalCells = (AliAODCaloCells*)fAODObjects->FindObject("emcalCells");
251 fPhosCells = (AliAODCaloCells*)fAODObjects->FindObject("phosCells");
252 fCaloClusters = (TClonesArray*)fAODObjects->FindObject("caloClusters");
253 fFmdClusters = (TClonesArray*)fAODObjects->FindObject("fmdClusters");
254 fPmdClusters = (TClonesArray*)fAODObjects->FindObject("pmdClusters");
257 //______________________________________________________________________________
258 void AliAODEvent::ResetStd(Int_t trkArrSize,
266 // deletes content of standard arrays and resets size
269 if (trkArrSize > fTracks->GetSize())
270 fTracks->Expand(trkArrSize);
273 if (vtxArrSize > fVertices->GetSize())
274 fVertices->Expand(vtxArrSize);
277 if (v0ArrSize > fV0s->GetSize())
278 fV0s->Expand(v0ArrSize);
281 if (jetSize > fJets->GetSize())
282 fJets->Expand(jetSize);
284 fCaloClusters->Delete();
285 if (caloClusSize > fCaloClusters->GetSize())
286 fCaloClusters->Expand(caloClusSize);
288 fFmdClusters->Delete();
289 if (fmdClusSize > fFmdClusters->GetSize())
290 fFmdClusters->Expand(fmdClusSize);
292 fPmdClusters->Delete();
293 if (pmdClusSize > fPmdClusters->GetSize())
294 fPmdClusters->Expand(pmdClusSize);
296 // Reset the tracklets
297 fTracklets->DeleteContainer();
298 fPhosCells->DeleteContainer();
299 fEmcalCells->DeleteContainer();
303 void AliAODEvent::ClearStd()
305 // clears the standard arrays
309 fTracklets ->DeleteContainer();
311 fEmcalCells ->DeleteContainer();
312 fPhosCells ->DeleteContainer();
313 fCaloClusters ->Clear();
314 fFmdClusters ->Clear();
315 fPmdClusters ->Clear();
318 //______________________________________________________________________________
319 Int_t AliAODEvent::GetMuonTracks(TRefArray *muonTracks) const
321 // fills the provided TRefArray with all found muon tracks
325 AliAODTrack *track = 0;
326 for (Int_t iTrack = 0; iTrack < GetNTracks(); iTrack++) {
327 if ((track = GetTrack(iTrack))->IsMuonTrack()) {
328 muonTracks->Add(track);
332 return muonTracks->GetEntriesFast();
336 void AliAODEvent::ReadFromTree(TTree *tree)
338 // connects aod event to tree
341 Printf("%s %d AliAODEvent::ReadFromTree() Zero Pointer to Tree \n",(char*)__FILE__,__LINE__);
345 if(!tree->GetTree())tree->LoadTree(0);
347 // Try to find AliAODEvent
348 AliAODEvent *aodEvent = 0;
349 aodEvent = (AliAODEvent*)tree->GetTree()->GetUserInfo()->FindObject("AliAODEvent");
351 // Check if already connected to tree
352 TList* connectedList = (TList*) (tree->GetUserInfo()->FindObject("AODObjectsConnectedToTree"));
354 // If connected use the connected list if objects
355 fAODObjects->Delete();
356 fAODObjects = connectedList;
361 // prevent a memory leak when reading back the TList
364 // create a new TList from the UserInfo TList...
365 // copy constructor does not work...
366 fAODObjects = (TList*)(aodEvent->GetList()->Clone());
367 fAODObjects->SetOwner(kFALSE);
368 if(fAODObjects->GetEntries()<kAODListN){
369 printf("%s %d AliAODEvent::ReadFromTree() TList contains less than the standard contents %d < %d \n",
370 (char*)__FILE__,__LINE__,fAODObjects->GetEntries(),kAODListN);
372 // set the branch addresses
373 TIter next(fAODObjects);
375 while((el=(TNamed*)next())){
376 TString bname(el->GetName());
377 // check if branch exists under this Name
378 TBranch *br = tree->GetBranch(bname.Data());
380 tree->SetBranchAddress(bname.Data(),fAODObjects->GetObjectRef(el));
383 br = tree->GetBranch(Form("%s.",bname.Data()));
385 tree->SetBranchAddress(Form("%s.",bname.Data()),fAODObjects->GetObjectRef(el));
388 printf("%s %d AliAODEvent::ReadFromTree() No Branch found with Name %s. \n",
389 (char*)__FILE__,__LINE__,bname.Data());
395 // when reading back we are not owner of the list
396 // must not delete it
397 fAODObjects->SetOwner(kFALSE);
398 fAODObjects->SetName("AODObjectsConnectedToTree");
399 // we are not owner of the list objects
400 // must not delete it
401 tree->GetUserInfo()->Add(fAODObjects);
404 // we can't get the list from the user data, create standard content
405 // and set it by hand
407 TIter next(fAODObjects);
409 while((el=(TNamed*)next())){
410 TString bname(el->GetName());
411 tree->SetBranchAddress(bname.Data(),fAODObjects->GetObjectRef(el));
414 // when reading back we are not owner of the list
415 // must not delete it
416 fAODObjects->SetOwner(kFALSE);
420 //______________________________________________________________________________
421 void AliAODEvent::Print(Option_t *) const
423 // Something meaningful should be implemented here.