]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliAODEvent.cxx
Add inheritance from AliVEvent. Fix standard naming scheme. Fix to be able to ReadFro...
[u/mrichter/AliRoot.git] / STEER / AliAODEvent.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2007, 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 //     AOD base class
20 //     Author: Markus Oldenburg, CERN
21 //-------------------------------------------------------------------------
22
23 #include <TTree.h>
24
25 #include "AliAODEvent.h"
26 #include "AliAODHeader.h"
27 #include "AliAODTrack.h"
28
29 ClassImp(AliAODEvent)
30
31 // definition of std AOD member names
32   const char* AliAODEvent::fAODListName[kAODListN] = {"header",
33                                                       "tracks",
34                                                       "vertices",
35                                                       "clusters",
36                                                       "jets",
37                                                       "tracklets"};
38 //______________________________________________________________________________
39 AliAODEvent::AliAODEvent() :
40   AliVEvent(),
41   fAODObjects(new TList()),
42   fHeader(0),
43   fTracks(0),
44   fVertices(0),
45   fClusters(0),
46   fJets(0),
47   fTracklets(0)
48 {
49   // default constructor
50 }
51
52 //______________________________________________________________________________
53 AliAODEvent::~AliAODEvent() 
54 {
55 // destructor
56     delete fAODObjects;
57 }
58
59 //______________________________________________________________________________
60 void AliAODEvent::AddObject(TObject* obj) 
61 {
62   // Add an object to the list of object.
63   // Please be aware that in order to increase performance you should
64   // refrain from using TObjArrays (if possible). Use TClonesArrays, instead.
65   
66   fAODObjects->AddLast(obj);
67 }
68
69 //______________________________________________________________________________
70 TObject *AliAODEvent::FindListObject(const char *objName)
71 {
72   // Return the pointer to the object with the given name.
73
74   return fAODObjects->FindObject(objName);
75 }
76
77 //______________________________________________________________________________
78 void AliAODEvent::CreateStdContent() 
79 {
80   // create the standard AOD content and set pointers
81
82   // create standard objects and add them to the TList of objects
83   AddObject(new AliAODHeader());
84   AddObject(new TClonesArray("AliAODTrack", 0));
85   AddObject(new TClonesArray("AliAODVertex", 0));
86   AddObject(new TClonesArray("AliAODCluster", 0));
87   AddObject(new TClonesArray("AliAODJet", 0));
88   AddObject(new AliAODTracklets());
89
90   // set names
91   SetStdNames();
92
93   // read back pointers
94   GetStdContent();
95
96   return;
97 }
98
99 //______________________________________________________________________________
100 void AliAODEvent::SetStdNames()
101 {
102   // introduce the standard naming
103
104   if(fAODObjects->GetEntries()==kAODListN){
105     for(int i = 0;i < fAODObjects->GetEntries();i++){
106       TObject *fObj = fAODObjects->At(i);
107       if(fObj->InheritsFrom("TNamed")){
108         ((TNamed*)fObj)->SetName(fAODListName[i]);
109       }
110       else if(fObj->InheritsFrom("TClonesArray")){
111         ((TClonesArray*)fObj)->SetName(fAODListName[i]);
112       }
113     }
114   }
115   else{
116     printf("%s:%d SetStdNames() Wrong number of Std Entries \n",(char*)__FILE__,__LINE__);
117   }
118
119
120 //______________________________________________________________________________
121 void AliAODEvent::GetStdContent()
122 {
123   // set pointers for standard content
124
125   fHeader    = (AliAODHeader*)fAODObjects->FindObject("header");
126   fTracks    = (TClonesArray*)fAODObjects->FindObject("tracks");
127   fVertices  = (TClonesArray*)fAODObjects->FindObject("vertices");
128   fClusters  = (TClonesArray*)fAODObjects->FindObject("clusters");
129   fJets      = (TClonesArray*)fAODObjects->FindObject("jets");
130   fTracklets = (AliAODTracklets*)fAODObjects->FindObject("tracklets");
131   
132   /* old implementation
133   fHeader    = (AliAODHeader*)fAODObjects->At(0);
134   fTracks    = (TClonesArray*)fAODObjects->At(1);
135   fVertices  = (TClonesArray*)fAODObjects->At(2);
136   fClusters  = (TClonesArray*)fAODObjects->At(3);
137   fJets      = (TClonesArray*)fAODObjects->At(4);
138   fTracklets = (AliAODTracklets*)fAODObjects->At(5);
139   */
140 }
141
142 //______________________________________________________________________________
143 void AliAODEvent::ResetStd(Int_t trkArrSize, Int_t vtxArrSize)
144 {
145   // deletes content of standard arrays and resets size
146   fTracks->Delete();
147   if (trkArrSize > fTracks->GetSize()) 
148     fTracks->Expand(trkArrSize);
149
150   fVertices->Delete();
151   if (vtxArrSize > fVertices->GetSize()) 
152     fVertices->Expand(vtxArrSize);
153 }
154
155 void AliAODEvent::ClearStd()
156 {
157   // clears the standard arrays
158     fTracks   ->Clear();
159     fVertices ->Clear();
160     fClusters ->Clear();
161     fJets     ->Clear();
162     fTracklets->DeleteContainer();
163 }
164
165 //______________________________________________________________________________
166 Int_t AliAODEvent::GetMuonTracks(TRefArray *muonTracks) const
167 {
168   // fills the provided TRefArray with all found muon tracks
169
170   muonTracks->Clear();
171
172   AliAODTrack *track = 0;
173   for (Int_t iTrack = 0; iTrack < GetNTracks(); iTrack++) {
174     if ((track = GetTrack(iTrack))->IsMuonTrack()) {
175       muonTracks->Add(track);
176     }
177   }
178   
179   return muonTracks->GetSize();
180 }
181
182
183 void AliAODEvent::ReadFromTree(TTree *tree)
184 {
185   // connects aod event to tree
186   
187   // load the TTree
188   tree->LoadTree(0);
189
190   fAODObjects = (TList*)((AliAODEvent*)tree->GetTree()->GetUserInfo()->FindObject("AliAODEvent"))->GetList(); 
191   TIter next(fAODObjects);
192   TNamed *el;
193   while((el=(TNamed*)next())){
194     TString bname(el->GetName());
195     tree->SetBranchAddress(bname.Data(),fAODObjects->GetObjectRef(el));
196   }
197   GetStdContent();
198 }
199
200 //______________________________________________________________________________
201 void AliAODEvent::Print(Option_t *) const
202 {
203   // Something meaningful should be implemented here.
204   
205   return;
206 }