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>
27 #include <TProcessID.h>
28 #include <TCollection.h>
29 #include "Riostream.h"
30 #include "AliAODEvent.h"
31 #include "AliAODHeader.h"
32 #include "AliAODTrack.h"
36 // definition of std AOD member names
37 const char* AliAODEvent::fAODListName[kAODListN] = {"header",
51 //______________________________________________________________________________
52 AliAODEvent::AliAODEvent() :
54 fAODObjects(new TList()),
70 // default constructor
73 //______________________________________________________________________________
74 AliAODEvent::AliAODEvent(const AliAODEvent& aod):
76 fAODObjects(new TList()),
77 fAODFolder(new TFolder()),
79 fHeader(new AliAODHeader(*aod.fHeader)),
80 fTracks(new TClonesArray(*aod.fTracks)),
81 fVertices(new TClonesArray(*aod.fVertices)),
82 fV0s(new TClonesArray(*aod.fV0s)),
83 fCascades(new TClonesArray(*aod.fCascades)),
84 fTracklets(new AliAODTracklets(*aod.fTracklets)),
85 fJets(new TClonesArray(*aod.fJets)),
86 fEmcalCells(new AliAODCaloCells(*aod.fEmcalCells)),
87 fPhosCells(new AliAODCaloCells(*aod.fPhosCells)),
88 fCaloClusters(new TClonesArray(*aod.fCaloClusters)),
89 fFmdClusters(new TClonesArray(*aod.fFmdClusters)),
90 fPmdClusters(new TClonesArray(*aod.fPmdClusters))
98 AddObject(fTracklets);
100 AddObject(fEmcalCells);
101 AddObject(fPhosCells);
102 AddObject(fCaloClusters);
103 AddObject(fFmdClusters);
104 AddObject(fPmdClusters);
105 fConnected = aod.fConnected;
109 //______________________________________________________________________________
110 AliAODEvent & AliAODEvent::operator=(const AliAODEvent& aod) {
112 // Assignment operator
114 if(&aod == this) return *this;
115 AliVEvent::operator=(aod);
117 fAODObjects = new TList();
118 fAODFolder = new TFolder();
119 fConnected = aod.fConnected;
120 fHeader = new AliAODHeader(*aod.fHeader);
121 fTracks = new TClonesArray(*aod.fTracks);
122 fVertices = new TClonesArray(*aod.fVertices);
123 fV0s = new TClonesArray(*aod.fV0s);
124 fCascades = new TClonesArray(*aod.fCascades);
125 fTracklets = new AliAODTracklets(*aod.fTracklets);
126 fJets = new TClonesArray(*aod.fJets);
127 fEmcalCells = new AliAODCaloCells(*aod.fEmcalCells);
128 fPhosCells = new AliAODCaloCells(*aod.fPhosCells);
129 fCaloClusters = new TClonesArray(*aod.fCaloClusters);
130 fFmdClusters = new TClonesArray(*aod.fFmdClusters);
131 fPmdClusters = new TClonesArray(*aod.fPmdClusters);
136 AddObject(fVertices);
138 AddObject(fCascades);
139 AddObject(fTracklets);
141 AddObject(fEmcalCells);
142 AddObject(fPhosCells);
143 AddObject(fCaloClusters);
144 AddObject(fFmdClusters);
145 AddObject(fPmdClusters);
151 //______________________________________________________________________________
152 AliAODEvent::~AliAODEvent()
155 if(fAODObjects&&!fConnected)
164 //______________________________________________________________________________
165 void AliAODEvent::AddObject(TObject* obj)
167 // Add an object to the list of objects.
168 // Please be aware that in order to increase performance you should
169 // refrain from using TObjArrays (if possible). Use TClonesArrays, instead.
171 fAODObjects->AddLast(obj);
174 //______________________________________________________________________________
175 void AliAODEvent::RemoveObject(TObject* obj)
177 // Removes an object from the list of objects.
179 fAODObjects->Remove(obj);
182 //______________________________________________________________________________
183 TObject *AliAODEvent::FindListObject(const char *objName)
185 // Return the pointer to the object with the given name.
187 return fAODObjects->FindObject(objName);
190 //______________________________________________________________________________
191 void AliAODEvent::CreateStdContent()
193 // create the standard AOD content and set pointers
195 // create standard objects and add them to the TList of objects
196 AddObject(new AliAODHeader());
197 AddObject(new TClonesArray("AliAODTrack", 0));
198 AddObject(new TClonesArray("AliAODVertex", 0));
199 AddObject(new TClonesArray("AliAODv0", 0));
200 AddObject(new TClonesArray("AliAODcascade", 0));
201 AddObject(new AliAODTracklets());
202 AddObject(new TClonesArray("AliAODJet", 0));
203 AddObject(new AliAODCaloCells());
204 AddObject(new AliAODCaloCells());
205 AddObject(new TClonesArray("AliAODCaloCluster", 0));
206 AddObject(new TClonesArray("AliAODFmdCluster", 0));
207 AddObject(new TClonesArray("AliAODPmdCluster", 0));
211 // read back pointers
217 void AliAODEvent::MakeEntriesReferencable()
219 // Make all entries referencable in a subsequent process
221 TIter next(fAODObjects);
223 while ((obj = next()))
225 if(obj->InheritsFrom("TCollection"))
227 AssignIDtoCollection((TCollection*)obj);
232 //______________________________________________________________________________
233 void AliAODEvent::SetStdNames()
235 // introduce the standard naming
237 if(fAODObjects->GetEntries()==kAODListN){
238 for(int i = 0;i < fAODObjects->GetEntries();i++){
239 TObject *fObj = fAODObjects->At(i);
240 if(fObj->InheritsFrom("TNamed")){
241 ((TNamed*)fObj)->SetName(fAODListName[i]);
243 else if(fObj->InheritsFrom("TClonesArray")){
244 ((TClonesArray*)fObj)->SetName(fAODListName[i]);
249 printf("%s:%d SetStdNames() Wrong number of Std Entries \n",(char*)__FILE__,__LINE__);
253 void AliAODEvent::CreateStdFolders()
255 // Create the standard folder structure
256 fAODFolder = gROOT->GetRootFolder()->AddFolder("AOD", "AOD");
257 if(fAODObjects->GetEntries()==kAODListN){
258 for(int i = 0;i < fAODObjects->GetEntries();i++){
259 TObject *fObj = fAODObjects->At(i);
260 if(fObj->InheritsFrom("TClonesArray")){
261 fAODFolder->AddFolder(fAODListName[i], fAODListName[i], (TCollection*) fObj);
263 fAODFolder->AddFolder(fAODListName[i], fAODListName[i], 0);
268 printf("%s:%d CreateStdFolders() Wrong number of Std Entries \n",(char*)__FILE__,__LINE__);
272 //______________________________________________________________________________
273 void AliAODEvent::GetStdContent()
275 // set pointers for standard content
277 fHeader = (AliAODHeader*)fAODObjects->FindObject("header");
278 fTracks = (TClonesArray*)fAODObjects->FindObject("tracks");
279 fVertices = (TClonesArray*)fAODObjects->FindObject("vertices");
280 fV0s = (TClonesArray*)fAODObjects->FindObject("v0s");
281 fCascades = (TClonesArray*)fAODObjects->FindObject("cascades");
282 fTracklets = (AliAODTracklets*)fAODObjects->FindObject("tracklets");
283 fJets = (TClonesArray*)fAODObjects->FindObject("jets");
284 fEmcalCells = (AliAODCaloCells*)fAODObjects->FindObject("emcalCells");
285 fPhosCells = (AliAODCaloCells*)fAODObjects->FindObject("phosCells");
286 fCaloClusters = (TClonesArray*)fAODObjects->FindObject("caloClusters");
287 fFmdClusters = (TClonesArray*)fAODObjects->FindObject("fmdClusters");
288 fPmdClusters = (TClonesArray*)fAODObjects->FindObject("pmdClusters");
291 //______________________________________________________________________________
292 void AliAODEvent::ResetStd(Int_t trkArrSize,
295 Int_t cascadeArrSize,
302 // deletes content of standard arrays and resets size
305 if (trkArrSize > fTracks->GetSize())
306 fTracks->Expand(trkArrSize);
309 if (vtxArrSize > fVertices->GetSize())
310 fVertices->Expand(vtxArrSize);
313 if (v0ArrSize > fV0s->GetSize())
314 fV0s->Expand(v0ArrSize);
317 if (cascadeArrSize > fCascades->GetSize())
318 fCascades->Expand(cascadeArrSize);
321 if (jetSize > fJets->GetSize())
322 fJets->Expand(jetSize);
324 fCaloClusters->Delete();
325 if (caloClusSize > fCaloClusters->GetSize())
326 fCaloClusters->Expand(caloClusSize);
328 fFmdClusters->Delete();
329 if (fmdClusSize > fFmdClusters->GetSize())
330 fFmdClusters->Expand(fmdClusSize);
332 fPmdClusters->Delete();
333 if (pmdClusSize > fPmdClusters->GetSize())
334 fPmdClusters->Expand(pmdClusSize);
336 // Reset the tracklets
337 fTracklets->DeleteContainer();
338 fPhosCells->DeleteContainer();
339 fEmcalCells->DeleteContainer();
343 void AliAODEvent::ClearStd()
345 // clears the standard arrays
346 fHeader ->RemoveQTheta();
348 fVertices ->Delete();
350 fCascades ->Delete();
351 fTracklets ->DeleteContainer();
353 fEmcalCells ->DeleteContainer();
354 fPhosCells ->DeleteContainer();
355 fCaloClusters ->Delete();
356 fFmdClusters ->Clear();
357 fPmdClusters ->Clear();
360 //_________________________________________________________________
361 Int_t AliAODEvent::GetPHOSClusters(TRefArray *clusters) const
363 // fills the provided TRefArray with all found phos clusters
367 AliAODCaloCluster *cl = 0;
368 Bool_t first = kTRUE;
369 for (Int_t i = 0; i < GetNCaloClusters() ; i++) {
370 if ( (cl = GetCaloCluster(i)) ) {
371 if (cl->IsPHOSCluster()){
373 new (clusters) TRefArray(TProcessID::GetProcessWithUID(cl));
377 //printf("IsPHOS cluster %d, E %2.3f Size: %d \n",i,cl->E(),clusters->GetEntriesFast());
381 return clusters->GetEntriesFast();
384 //_________________________________________________________________
385 Int_t AliAODEvent::GetEMCALClusters(TRefArray *clusters) const
387 // fills the provided TRefArray with all found emcal clusters
390 AliAODCaloCluster *cl = 0;
391 Bool_t first = kTRUE;
392 for (Int_t i = 0; i < GetNCaloClusters(); i++) {
393 if ( (cl = GetCaloCluster(i)) ) {
394 if (cl->IsEMCALCluster()){
396 new (clusters) TRefArray(TProcessID::GetProcessWithUID(cl));
400 //printf("IsEMCal cluster %d, E %2.3f Size: %d \n",i,cl->E(),clusters->GetEntriesFast());
404 return clusters->GetEntriesFast();
408 //______________________________________________________________________________
409 Int_t AliAODEvent::GetMuonTracks(TRefArray *muonTracks) const
411 // fills the provided TRefArray with all found muon tracks
415 AliAODTrack *track = 0;
416 for (Int_t iTrack = 0; iTrack < GetNTracks(); iTrack++) {
417 if ((track = GetTrack(iTrack))->IsMuonTrack()) {
418 muonTracks->Add(track);
422 return muonTracks->GetEntriesFast();
426 void AliAODEvent::ReadFromTree(TTree *tree, Option_t* opt /*= ""*/)
428 // Connects aod event to tree
431 Printf("%s %d AliAODEvent::ReadFromTree() Zero Pointer to Tree \n",(char*)__FILE__,__LINE__);
435 if(!tree->GetTree())tree->LoadTree(0);
437 // Try to find AliAODEvent
438 AliAODEvent *aodEvent = 0;
439 aodEvent = (AliAODEvent*)tree->GetTree()->GetUserInfo()->FindObject("AliAODEvent");
441 // Check if already connected to tree
442 TList* connectedList = (TList*) (tree->GetUserInfo()->FindObject("AODObjectsConnectedToTree"));
443 if (connectedList && (strcmp(opt, "reconnect"))) {
444 // If connected use the connected list of objects
445 fAODObjects->Delete();
446 fAODObjects = connectedList;
452 // prevent a memory leak when reading back the TList
453 if (!(strcmp(opt, "reconnect"))) fAODObjects->Delete();
455 // create a new TList from the UserInfo TList...
456 // copy constructor does not work...
457 fAODObjects = (TList*)(aodEvent->GetList()->Clone());
458 fAODObjects->SetOwner(kTRUE);
459 if(fAODObjects->GetEntries()<kAODListN){
460 printf("%s %d AliAODEvent::ReadFromTree() TList contains less than the standard contents %d < %d \n",
461 (char*)__FILE__,__LINE__,fAODObjects->GetEntries(),kAODListN);
464 // Let's find out whether we have friends
465 TList* friendL = tree->GetTree()->GetListOfFriends();
470 while ((fe = (TFriendElement*)next())){
471 aodEvent = (AliAODEvent*)(fe->GetTree()->GetUserInfo()->FindObject("AliAODEvent"));
473 printf("No UserInfo on tree \n");
476 TList* objL = (TList*)(aodEvent->GetList()->Clone());
477 printf("Get list of object from tree %d !!\n", objL->GetEntries());
478 TIter nextobject(objL);
480 while((obj = nextobject()))
482 printf("Adding object from friend %s !\n", obj->GetName());
483 fAODObjects->Add(obj);
484 } // object "branch" loop
490 // set the branch addresses
491 TIter next(fAODObjects);
493 while((el=(TNamed*)next())){
494 TString bname(el->GetName());
495 // check if branch exists under this Name
496 TBranch *br = tree->GetTree()->GetBranch(bname.Data());
498 tree->SetBranchAddress(bname.Data(),fAODObjects->GetObjectRef(el));
500 br = tree->GetBranch(Form("%s.",bname.Data()));
502 tree->SetBranchAddress(Form("%s.",bname.Data()),fAODObjects->GetObjectRef(el));
505 printf("%s %d AliAODEvent::ReadFromTree() No Branch found with Name %s. \n",
506 (char*)__FILE__,__LINE__,bname.Data());
511 // when reading back we are not owner of the list
512 // must not delete it
513 fAODObjects->SetOwner(kTRUE);
514 fAODObjects->SetName("AODObjectsConnectedToTree");
515 // we are not owner of the list objects
516 // must not delete it
517 tree->GetUserInfo()->Add(fAODObjects);
521 // we can't get the list from the user data, create standard content
522 // and set it by hand
524 TIter next(fAODObjects);
526 while((el=(TNamed*)next())){
527 TString bname(el->GetName());
528 tree->SetBranchAddress(bname.Data(),fAODObjects->GetObjectRef(el));
531 // when reading back we are not owner of the list
532 // must not delete it
533 fAODObjects->SetOwner(kTRUE);
537 //______________________________________________________________________________
538 void AliAODEvent::Print(Option_t *) const
540 // Print the names of the all branches
541 TIter next(fAODObjects);
543 Printf(">>>>> AOD Content <<<<<");
544 while((el=(TNamed*)next())){
545 Printf(">> %s ",el->GetName());
547 Printf(">>>>> <<<<<");
552 void AliAODEvent::AssignIDtoCollection(TCollection* col)
554 // Static method which assigns a ID to each object in a collection
555 // In this way the objects are marked as referenced and written with
556 // an ID. This has the advantage that TRefs to this objects can be
557 // written by a subsequent process.
560 while ((obj = next()))
561 TProcessID::AssignID(obj);