Line data Source code
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 <TROOT.h>
24 : #include <TTree.h>
25 : #include <TFolder.h>
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"
33 : #include "AliAODDimuon.h"
34 : #include "AliAODTrdTrack.h"
35 : #include "event.h"
36 :
37 170 : ClassImp(AliAODEvent)
38 :
39 : // definition of std AOD member names
40 : const char* AliAODEvent::fAODListName[kAODListN] = {"header",
41 : "tracks",
42 : "vertices",
43 : "v0s",
44 : "cascades",
45 : "tracklets",
46 : "jets",
47 : "emcalCells",
48 : "phosCells",
49 : "caloClusters",
50 : "emcalTrigger",
51 : "phosTrigger",
52 : "fmdClusters",
53 : "pmdClusters",
54 : "hmpidRings",
55 : "dimuons",
56 : "AliAODTZERO",
57 : "AliAODVZERO",
58 : "AliAODZDC",
59 : "AliAODAD",
60 : "AliTOFHeader",
61 : "trdTracks"
62 :
63 : };
64 : //______________________________________________________________________________
65 : AliAODEvent::AliAODEvent() :
66 2 : AliVEvent(),
67 2 : fAODObjects(0),
68 2 : fAODFolder(0),
69 2 : fConnected(kFALSE),
70 2 : fTracksConnected(kFALSE),
71 2 : fHeader(0),
72 2 : fTracks(0),
73 2 : fVertices(0),
74 2 : fV0s(0),
75 2 : fCascades(0),
76 2 : fTracklets(0),
77 2 : fJets(0),
78 2 : fEmcalCells(0),
79 2 : fPhosCells(0),
80 2 : fCaloClusters(0),
81 2 : fEMCALTrigger(0),
82 2 : fPHOSTrigger(0),
83 2 : fFmdClusters(0),
84 2 : fPmdClusters(0),
85 2 : fHMPIDrings(0),
86 2 : fDimuons(0),
87 2 : fAODTZERO(0),
88 2 : fAODVZERO(0),
89 2 : fAODZDC(0),
90 2 : fAODAD(0),
91 2 : fTOFHeader(0),
92 2 : fTrdTracks(0)
93 10 : {
94 : // default constructor
95 10 : if (TClass::IsCallingNew() != TClass::kDummyNew) fAODObjects = new TList();
96 4 : }
97 :
98 : //______________________________________________________________________________
99 : AliAODEvent::AliAODEvent(const AliAODEvent& aod):
100 0 : AliVEvent(aod),
101 0 : fAODObjects(new TList()),
102 0 : fAODFolder(0),
103 0 : fConnected(kFALSE),
104 0 : fTracksConnected(kFALSE),
105 : // fHeader(new AliAODHeader(*aod.fHeader)),
106 0 : fHeader(0),
107 0 : fTracks(new TClonesArray(*aod.fTracks)),
108 0 : fVertices(new TClonesArray(*aod.fVertices)),
109 0 : fV0s(new TClonesArray(*aod.fV0s)),
110 0 : fCascades(new TClonesArray(*aod.fCascades)),
111 0 : fTracklets(new AliAODTracklets(*aod.fTracklets)),
112 0 : fJets(new TClonesArray(*aod.fJets)),
113 0 : fEmcalCells(new AliAODCaloCells(*aod.fEmcalCells)),
114 0 : fPhosCells(new AliAODCaloCells(*aod.fPhosCells)),
115 0 : fCaloClusters(new TClonesArray(*aod.fCaloClusters)),
116 0 : fEMCALTrigger(new AliAODCaloTrigger(*aod.fEMCALTrigger)),
117 0 : fPHOSTrigger(new AliAODCaloTrigger(*aod.fPHOSTrigger)),
118 0 : fFmdClusters(new TClonesArray(*aod.fFmdClusters)),
119 0 : fPmdClusters(new TClonesArray(*aod.fPmdClusters)),
120 0 : fHMPIDrings(new TClonesArray(*aod.fHMPIDrings)),
121 0 : fDimuons(new TClonesArray(*aod.fDimuons)),
122 0 : fAODTZERO(new AliAODTZERO(*aod.fAODTZERO)),
123 0 : fAODVZERO(new AliAODVZERO(*aod.fAODVZERO)),
124 0 : fAODZDC(new AliAODZDC(*aod.fAODZDC)),
125 0 : fAODAD(new AliAODAD(*aod.fAODAD)),
126 0 : fTOFHeader(new AliTOFHeader(*aod.fTOFHeader)),
127 0 : fTrdTracks(new TClonesArray(*aod.fTrdTracks))
128 0 : {
129 : // Copy constructor
130 0 : AddHeader(fHeader);
131 0 : AddObject(fHeader);
132 0 : AddObject(fTracks);
133 0 : AddObject(fVertices);
134 0 : AddObject(fV0s);
135 0 : AddObject(fCascades);
136 0 : AddObject(fTracklets);
137 0 : AddObject(fJets);
138 0 : AddObject(fEmcalCells);
139 0 : AddObject(fPhosCells);
140 0 : AddObject(fCaloClusters);
141 0 : AddObject(fEMCALTrigger);
142 0 : AddObject(fPHOSTrigger);
143 0 : AddObject(fFmdClusters);
144 0 : AddObject(fPmdClusters);
145 0 : AddObject(fHMPIDrings);
146 0 : AddObject(fDimuons);
147 0 : AddObject(fAODTZERO);
148 0 : AddObject(fAODVZERO);
149 0 : AddObject(fAODZDC);
150 0 : AddObject(fAODAD);
151 0 : AddObject(fTOFHeader);
152 0 : AddObject(fTrdTracks);
153 0 : fConnected = aod.fConnected;
154 0 : ConnectTracks();
155 0 : GetStdContent();
156 0 : CreateStdFolders();
157 0 : }
158 :
159 : //______________________________________________________________________________
160 : AliAODEvent & AliAODEvent::operator=(const AliAODEvent& aod) {
161 :
162 : // Assignment operator
163 :
164 0 : if(&aod == this) return *this;
165 0 : AliVEvent::operator=(aod);
166 :
167 : // This assumes that the list is already created
168 : // and that the virtual void Copy(Tobject&) function
169 : // is correctly implemented in the derived class
170 : // otherwise only TObject::Copy() will be used
171 :
172 0 : if((fAODObjects->GetSize()==0)&&(aod.fAODObjects->GetSize()>=kAODListN)){
173 : // We cover the case that we do not yet have the
174 : // standard content but the source has it
175 0 : CreateStdContent();
176 0 : }
177 :
178 : // Here we have the standard content without user additions, but the content is
179 : // not matching the aod source.
180 :
181 : // Iterate the list of source objects
182 0 : TIter next(aod.GetList());
183 : TObject *its = 0;
184 0 : TString name;
185 0 : while ((its = next())) {
186 0 : name = its->GetName();
187 : // Check if we have this object type in out list
188 0 : TObject *mine = fAODObjects->FindObject(name);
189 0 : if(!mine) {
190 : // We have to create the same type of object.
191 0 : TClass* pClass=TClass::GetClass(its->ClassName());
192 0 : if (!pClass) {
193 0 : AliWarning(Form("Can not find class description for entry %s (%s)\n",
194 : its->ClassName(), name.Data()));
195 0 : continue;
196 : }
197 0 : mine=(TObject*)pClass->New();
198 0 : if(!mine){
199 : // not in this: can be added to list
200 0 : AliWarning(Form("%s:%d Could not find %s for copying \n",
201 : (char*)__FILE__,__LINE__,name.Data()));
202 0 : continue;
203 : }
204 0 : if(mine->InheritsFrom("TNamed")) {
205 0 : ((TNamed*)mine)->SetName(name);
206 0 : } else if(mine->InheritsFrom("TCollection")){
207 0 : if(mine->InheritsFrom("TClonesArray")) {
208 0 : TClonesArray *itscl = dynamic_cast<TClonesArray*>(its);
209 0 : if (!itscl) {
210 0 : AliWarning(Form("Class description for entry %s (%s) not TClonesArray\n",
211 : its->ClassName(), name.Data()));
212 0 : continue;
213 :
214 : }
215 0 : dynamic_cast<TClonesArray*>(mine)->SetClass(itscl->GetClass(), itscl->GetSize());
216 0 : }
217 0 : dynamic_cast<TCollection*>(mine)->SetName(name);
218 : }
219 0 : AliDebug(1, Form("adding object %s of type %s", mine->GetName(), mine->ClassName()));
220 0 : AddObject(mine);
221 0 : }
222 : // Now we have an object of the same type and name, but different content.
223 0 : if(!its->InheritsFrom("TCollection")){
224 : // simple objects (do they have a Copy method that calls operator= ?)
225 0 : its->Copy(*mine);
226 0 : } else if (its->InheritsFrom("TClonesArray")) {
227 : // Create or expand the tclonesarray pointers
228 : // so we can directly copy to the object
229 0 : TClonesArray *its_tca = (TClonesArray*)its;
230 0 : TClonesArray *mine_tca = (TClonesArray*)mine;
231 : // this leaves the capacity of the TClonesArray the same
232 : // except for a factor of 2 increase when size > capacity
233 : // does not release any memory occupied by the tca
234 0 : Int_t its_entries = its_tca->GetEntriesFast();
235 0 : mine_tca->ExpandCreate(its_entries);
236 0 : for(int i=0; i<its_entries; i++){
237 : // copy
238 0 : TObject *mine_tca_obj = mine_tca->At(i);
239 0 : TObject *its_tca_obj = its_tca->At(i);
240 : // no need to delete first
241 : // pointers within the class should be handled by Copy()...
242 : // Can there be Empty slots?
243 0 : its_tca_obj->Copy(*mine_tca_obj);
244 : }
245 0 : } else {
246 0 : AliWarning(Form("%s:%d cannot copy TCollection \n",
247 : (char*)__FILE__,__LINE__));
248 : }
249 0 : }
250 0 : fConnected = aod.fConnected;
251 0 : fTracksConnected = kFALSE;
252 0 : ConnectTracks();
253 : return *this;
254 0 : }
255 :
256 :
257 : //______________________________________________________________________________
258 : AliAODEvent::~AliAODEvent()
259 12 : {
260 : // destructor
261 4 : delete fAODFolder;
262 2 : fAODFolder = 0;
263 2 : if(!fConnected) {
264 : // fAODObjects->Delete("slow");
265 4 : delete fAODObjects;
266 : }
267 6 : }
268 :
269 : //______________________________________________________________________________
270 : void AliAODEvent::AddObject(TObject* obj)
271 : {
272 : // Add an object to the list of objects.
273 : // Please be aware that in order to increase performance you should
274 : // refrain from using TObjArrays (if possible). Use TClonesArrays, instead.
275 :
276 : // if ( !fAODObjects ) {
277 : // fAODObjects = new TList();
278 : // fAODObjects->SetOwner();
279 : // }
280 92 : if ( !fAODObjects->FindObject(obj) )
281 : {
282 46 : fAODObjects->AddLast(obj);
283 46 : }
284 46 : }
285 :
286 : //______________________________________________________________________________
287 : Int_t AliAODEvent::AddTrack(const AliAODTrack* trk)
288 : {
289 : // Add new AOD track. Make sure to set the event if needed.
290 0 : AliAODTrack *track = new((*fTracks)[fTracks->GetEntriesFast()]) AliAODTrack(*trk);
291 0 : track->SetAODEvent(this);
292 0 : return fTracks->GetEntriesFast()-1;
293 0 : }
294 :
295 : //______________________________________________________________________________
296 : void AliAODEvent::RemoveObject(TObject* obj)
297 : {
298 : // Removes an object from the list of objects.
299 :
300 0 : fAODObjects->Remove(obj);
301 0 : }
302 :
303 : //______________________________________________________________________________
304 : TObject *AliAODEvent::FindListObject(const char *objName) const
305 : {
306 : // Return the pointer to the object with the given name.
307 :
308 48 : return fAODObjects->FindObject(objName);
309 : }
310 :
311 : //______________________________________________________________________________
312 : void AliAODEvent::CreateStdContent()
313 : {
314 : // create the standard AOD content and set pointers
315 :
316 : // create standard objects and add them to the TList of objects
317 6 : AddObject(new AliAODHeader());
318 4 : AddObject(new TClonesArray("AliAODTrack", 0));
319 4 : AddObject(new TClonesArray("AliAODVertex", 0));
320 4 : AddObject(new TClonesArray("AliAODv0", 0));
321 4 : AddObject(new TClonesArray("AliAODcascade", 0));
322 4 : AddObject(new AliAODTracklets());
323 4 : AddObject(new TClonesArray("AliAODJet", 0));
324 4 : AddObject(new AliAODCaloCells());
325 4 : AddObject(new AliAODCaloCells());
326 4 : AddObject(new TClonesArray("AliAODCaloCluster", 0));
327 4 : AddObject(new AliAODCaloTrigger()); // EMCAL
328 4 : AddObject(new AliAODCaloTrigger()); // PHOS
329 4 : AddObject(new TClonesArray("AliAODFmdCluster", 0));
330 4 : AddObject(new TClonesArray("AliAODPmdCluster", 0));
331 4 : AddObject(new TClonesArray("AliAODHMPIDrings", 0));
332 4 : AddObject(new TClonesArray("AliAODDimuon", 0));
333 4 : AddObject(new AliAODTZERO());
334 4 : AddObject(new AliAODVZERO());
335 4 : AddObject(new AliAODZDC());
336 4 : AddObject(new AliAODAD());
337 4 : AddObject(new AliTOFHeader());
338 4 : AddObject(new TClonesArray("AliAODTrdTrack", 0));
339 : // set names
340 2 : SetStdNames();
341 :
342 : // read back pointers
343 2 : GetStdContent();
344 2 : CreateStdFolders();
345 2 : return;
346 0 : }
347 :
348 : void AliAODEvent::MakeEntriesReferencable()
349 : {
350 : // Make all entries referencable in a subsequent process
351 : //
352 16 : TIter next(fAODObjects);
353 : TObject* obj;
354 392 : while ((obj = next()))
355 : {
356 368 : if(obj->InheritsFrom("TCollection"))
357 : {
358 92 : AssignIDtoCollection((TCollection*)obj);
359 : }
360 : }
361 8 : }
362 :
363 : //______________________________________________________________________________
364 : void AliAODEvent::SetStdNames()
365 : {
366 : // introduce the standard naming
367 :
368 4 : if(fAODObjects->GetEntries()==kAODListN){
369 92 : for(int i = 0;i < fAODObjects->GetEntries();i++){
370 44 : TObject *fObj = fAODObjects->At(i);
371 44 : if(fObj->InheritsFrom("TNamed")){
372 12 : ((TNamed*)fObj)->SetName(fAODListName[i]);
373 12 : }
374 32 : else if(fObj->InheritsFrom("TClonesArray")){
375 22 : ((TClonesArray*)fObj)->SetName(fAODListName[i]);
376 22 : }
377 : }
378 2 : }
379 : else{
380 0 : printf("%s:%d SetStdNames() Wrong number of Std Entries \n",(char*)__FILE__,__LINE__);
381 : }
382 2 : }
383 :
384 : //______________________________________________________________________________
385 : void AliAODEvent::CreateStdFolders()
386 : {
387 : // Create the standard folder structure
388 4 : if(fAODFolder)delete fAODFolder;
389 2 : fAODFolder = gROOT->GetRootFolder()->AddFolder("AOD", "AOD");
390 2 : if(fAODObjects->GetEntries()==kAODListN){
391 92 : for(int i = 0;i < fAODObjects->GetEntries();i++){
392 44 : TObject *fObj = fAODObjects->At(i);
393 88 : if(fObj->InheritsFrom("TClonesArray")){
394 66 : fAODFolder->AddFolder(fAODListName[i], fAODListName[i], (TCollection*) fObj);
395 22 : } else {
396 22 : fAODFolder->AddFolder(fAODListName[i], fAODListName[i], 0);
397 : }
398 : }
399 2 : }
400 : else{
401 0 : printf("%s:%d CreateStdFolders() Wrong number of Std Entries \n",(char*)__FILE__,__LINE__);
402 : }
403 2 : }
404 :
405 : //______________________________________________________________________________
406 : void AliAODEvent::GetStdContent()
407 : {
408 : // set pointers for standard content
409 :
410 4 : fHeader = (AliVAODHeader*)fAODObjects->FindObject("header");
411 2 : fTracks = (TClonesArray*)fAODObjects->FindObject("tracks");
412 2 : fVertices = (TClonesArray*)fAODObjects->FindObject("vertices");
413 2 : fV0s = (TClonesArray*)fAODObjects->FindObject("v0s");
414 2 : fCascades = (TClonesArray*)fAODObjects->FindObject("cascades");
415 2 : fTracklets = (AliAODTracklets*)fAODObjects->FindObject("tracklets");
416 2 : fJets = (TClonesArray*)fAODObjects->FindObject("jets");
417 2 : fEmcalCells = (AliAODCaloCells*)fAODObjects->FindObject("emcalCells");
418 2 : fPhosCells = (AliAODCaloCells*)fAODObjects->FindObject("phosCells");
419 2 : fCaloClusters = (TClonesArray*)fAODObjects->FindObject("caloClusters");
420 2 : fEMCALTrigger = (AliAODCaloTrigger*)fAODObjects->FindObject("emcalTrigger");
421 2 : fPHOSTrigger = (AliAODCaloTrigger*)fAODObjects->FindObject("phosTrigger");
422 2 : fFmdClusters = (TClonesArray*)fAODObjects->FindObject("fmdClusters");
423 2 : fPmdClusters = (TClonesArray*)fAODObjects->FindObject("pmdClusters");
424 2 : fHMPIDrings = (TClonesArray*)fAODObjects->FindObject("hmpidRings");
425 2 : fDimuons = (TClonesArray*)fAODObjects->FindObject("dimuons");
426 2 : fAODTZERO = (AliAODTZERO*)fAODObjects->FindObject("AliAODTZERO");
427 2 : fAODVZERO = (AliAODVZERO*)fAODObjects->FindObject("AliAODVZERO");
428 2 : fAODZDC = (AliAODZDC*)fAODObjects->FindObject("AliAODZDC");
429 2 : fAODAD = (AliAODAD*)fAODObjects->FindObject("AliAODAD");
430 2 : fTOFHeader = (AliTOFHeader*)fAODObjects->FindObject("AliTOFHeader");
431 2 : fTrdTracks = (TClonesArray*)fAODObjects->FindObject("trdTracks");
432 2 : }
433 :
434 : //______________________________________________________________________________
435 : void AliAODEvent::ResetStd(Int_t trkArrSize,
436 : Int_t vtxArrSize,
437 : Int_t v0ArrSize,
438 : Int_t cascadeArrSize,
439 : Int_t jetSize,
440 : Int_t caloClusSize,
441 : Int_t fmdClusSize,
442 : Int_t pmdClusSize,
443 : Int_t hmpidRingsSize,
444 : Int_t dimuonArrSize,
445 : Int_t nTrdTracks
446 : )
447 : {
448 : // deletes content of standard arrays and resets size
449 32 : fTracksConnected = kFALSE;
450 16 : if (fTracks) {
451 16 : fTracks->Delete();
452 16 : if (trkArrSize > fTracks->GetSize())
453 2 : fTracks->Expand(trkArrSize);
454 : }
455 16 : if (fVertices) {
456 16 : fVertices->Delete();
457 16 : if (vtxArrSize > fVertices->GetSize())
458 0 : fVertices->Expand(vtxArrSize);
459 : }
460 16 : if (fV0s) {
461 16 : fV0s->Delete();
462 16 : if (v0ArrSize > fV0s->GetSize())
463 2 : fV0s->Expand(v0ArrSize);
464 : }
465 16 : if (fCascades) {
466 16 : fCascades->Delete();
467 16 : if (cascadeArrSize > fCascades->GetSize())
468 0 : fCascades->Expand(cascadeArrSize);
469 : }
470 16 : if (fJets) {
471 16 : fJets->Delete();
472 16 : if (jetSize > fJets->GetSize())
473 0 : fJets->Expand(jetSize);
474 : }
475 16 : if (fCaloClusters) {
476 16 : fCaloClusters->Delete();
477 16 : if (caloClusSize > fCaloClusters->GetSize())
478 1 : fCaloClusters->Expand(caloClusSize);
479 : }
480 16 : if (fFmdClusters) {
481 16 : fFmdClusters->Delete();
482 16 : if (fmdClusSize > fFmdClusters->GetSize())
483 0 : fFmdClusters->Expand(fmdClusSize);
484 : }
485 16 : if (fPmdClusters) {
486 16 : fPmdClusters->Delete();
487 16 : if (pmdClusSize > fPmdClusters->GetSize())
488 4 : fPmdClusters->Expand(pmdClusSize);
489 : }
490 16 : if (fHMPIDrings) {
491 16 : fHMPIDrings->Delete();
492 16 : if (hmpidRingsSize > fHMPIDrings->GetSize())
493 0 : fHMPIDrings->Expand(hmpidRingsSize);
494 : }
495 16 : if (fDimuons) {
496 16 : fDimuons->Delete();
497 16 : if (dimuonArrSize > fDimuons->GetSize())
498 0 : fDimuons->Expand(dimuonArrSize);
499 : }
500 16 : if (fTrdTracks) {
501 : // no pointers in there, so cheaper Clear suffices
502 : // fTrdTracks->Clear("C");
503 : // Not quite: AliAODTrdTrack has a clones array of tracklets inside
504 16 : fTrdTracks->Delete();
505 16 : if (nTrdTracks > fTrdTracks->GetSize())
506 0 : fTrdTracks->Expand(nTrdTracks);
507 : }
508 :
509 16 : if (fTracklets)
510 16 : fTracklets->DeleteContainer();
511 16 : if (fPhosCells)
512 16 : fPhosCells->DeleteContainer();
513 16 : if (fEmcalCells)
514 16 : fEmcalCells->DeleteContainer();
515 :
516 16 : if (fEMCALTrigger)
517 16 : fEMCALTrigger->DeAllocate();
518 16 : if (fPHOSTrigger)
519 16 : fPHOSTrigger->DeAllocate();
520 :
521 16 : }
522 :
523 : //______________________________________________________________________________
524 : void AliAODEvent::ClearStd()
525 : {
526 : // clears the standard arrays
527 0 : if (fHeader){
528 : // FIXME: this if-else patch was introduced by Michele Floris on 17/03/14 to test nano AOD. To be removed.
529 0 : if(fHeader->InheritsFrom("AliAODHeader")){
530 0 : fHeader ->Clear();
531 0 : }
532 : else {
533 : AliVHeader * head = 0;
534 0 : head = dynamic_cast<AliVHeader*>((TObject*)fHeader);
535 0 : if(head) head->Clear();
536 : }
537 : }
538 0 : fTracksConnected = kFALSE;
539 0 : if (fTracks)
540 0 : fTracks ->Delete();
541 0 : if (fVertices)
542 0 : fVertices ->Delete();
543 0 : if (fV0s)
544 0 : fV0s ->Delete();
545 0 : if (fCascades)
546 0 : fCascades ->Delete();
547 0 : if (fTracklets)
548 0 : fTracklets ->DeleteContainer();
549 0 : if (fJets)
550 0 : fJets ->Delete();
551 0 : if (fEmcalCells)
552 0 : fEmcalCells ->DeleteContainer();
553 0 : if (fPhosCells)
554 0 : fPhosCells ->DeleteContainer();
555 0 : if (fCaloClusters)
556 0 : fCaloClusters ->Delete();
557 0 : if (fFmdClusters)
558 0 : fFmdClusters ->Clear();
559 0 : if (fPmdClusters)
560 0 : fPmdClusters ->Clear();
561 0 : if (fHMPIDrings)
562 0 : fHMPIDrings ->Clear();
563 0 : if (fDimuons)
564 0 : fDimuons ->Clear("C");
565 0 : if (fTrdTracks)
566 0 : fTrdTracks ->Clear();
567 :
568 0 : if (fEMCALTrigger)
569 0 : fEMCALTrigger->DeAllocate();
570 0 : if (fPHOSTrigger)
571 0 : fPHOSTrigger->DeAllocate();
572 0 : }
573 :
574 : //_________________________________________________________________
575 : Int_t AliAODEvent::GetPHOSClusters(TRefArray *clusters) const
576 : {
577 : // fills the provided TRefArray with all found phos clusters
578 :
579 16 : clusters->Clear();
580 :
581 : AliAODCaloCluster *cl = 0;
582 : Bool_t first = kTRUE;
583 102 : for (Int_t i = 0; i < GetNumberOfCaloClusters() ; i++) {
584 43 : if ( (cl = GetCaloCluster(i)) ) {
585 43 : if (cl->IsPHOS()){
586 10 : if(first) {
587 16 : new (clusters) TRefArray(TProcessID::GetProcessWithUID(cl));
588 : first=kFALSE;
589 8 : }
590 10 : clusters->Add(cl);
591 : //printf("IsPHOS cluster %d, E %2.3f Size: %d \n",i,cl->E(),clusters->GetEntriesFast());
592 10 : }
593 : }
594 : }
595 8 : return clusters->GetEntriesFast();
596 0 : }
597 :
598 : //_________________________________________________________________
599 : Int_t AliAODEvent::GetEMCALClusters(TRefArray *clusters) const
600 : {
601 : // fills the provided TRefArray with all found emcal clusters
602 :
603 16 : clusters->Clear();
604 : AliAODCaloCluster *cl = 0;
605 : Bool_t first = kTRUE;
606 102 : for (Int_t i = 0; i < GetNumberOfCaloClusters(); i++) {
607 43 : if ( (cl = GetCaloCluster(i)) ) {
608 43 : if (cl->IsEMCAL()){
609 33 : if(first) {
610 16 : new (clusters) TRefArray(TProcessID::GetProcessWithUID(cl));
611 : first=kFALSE;
612 8 : }
613 33 : clusters->Add(cl);
614 : //printf("IsEMCal cluster %d, E %2.3f Size: %d \n",i,cl->E(),clusters->GetEntriesFast());
615 33 : }
616 : }
617 : }
618 8 : return clusters->GetEntriesFast();
619 0 : }
620 :
621 :
622 : //______________________________________________________________________________
623 : Int_t AliAODEvent::GetMuonTracks(TRefArray *muonTracks) const
624 : {
625 : // fills the provided TRefArray with all found muon tracks
626 :
627 0 : muonTracks->Clear();
628 :
629 : AliAODTrack *track = 0;
630 0 : for (Int_t iTrack = 0; iTrack < GetNumberOfTracks(); iTrack++) {
631 0 : track = dynamic_cast<AliAODTrack*>(GetTrack(iTrack));
632 0 : if(!track) AliFatal("Not a standard AOD");
633 0 : if (track->IsMuonTrack()) {
634 0 : muonTracks->Add(track);
635 0 : }
636 : }
637 :
638 0 : return muonTracks->GetEntriesFast();
639 : }
640 :
641 :
642 : //______________________________________________________________________________
643 : Int_t AliAODEvent::GetNumberOfMuonTracks() const
644 : {
645 : // get number of muon tracks
646 0 : int ntr = GetNumberOfTracks();
647 0 : if (!ntr) return 0;
648 : Int_t nMuonTracks=0;
649 0 : if(!dynamic_cast<AliAODTrack*>(GetTrack(0))) {
650 0 : AliError("Not a standard AOD");
651 0 : return 0;
652 : }
653 :
654 0 : for (Int_t iTrack=ntr; iTrack--;) {
655 0 : if (((AliAODTrack*)GetTrack(iTrack))->IsMuonTrack()) {
656 0 : nMuonTracks++;
657 0 : }
658 : }
659 :
660 0 : return nMuonTracks;
661 0 : }
662 :
663 : //______________________________________________________________________________
664 : Int_t AliAODEvent::GetMuonGlobalTracks(TRefArray *muonGlobalTracks) const // AU
665 : {
666 : // fills the provided TRefArray with all found muon global tracks
667 :
668 0 : muonGlobalTracks->Clear();
669 :
670 : AliAODTrack *track = 0;
671 0 : for (Int_t iTrack = 0; iTrack < GetNumberOfTracks(); iTrack++) {
672 0 : track = dynamic_cast<AliAODTrack*>(GetTrack(iTrack));
673 0 : if(!track) AliFatal("Not a standard AOD");
674 0 : if (track->IsMuonGlobalTrack()) {
675 0 : muonGlobalTracks->Add(track);
676 0 : }
677 : }
678 :
679 0 : return muonGlobalTracks->GetEntriesFast();
680 : }
681 :
682 :
683 : //______________________________________________________________________________
684 : Int_t AliAODEvent::GetNumberOfMuonGlobalTracks() const // AU
685 : {
686 : // get number of muon global tracks
687 0 : int ntr = GetNumberOfTracks();
688 0 : if (!ntr) return 0;
689 : Int_t nMuonGlobalTracks=0;
690 0 : if(!dynamic_cast<AliAODTrack*>(GetTrack(0))) {
691 0 : AliError("Not a standard AOD");
692 0 : return 0;
693 : }
694 0 : for (Int_t iTrack=ntr; iTrack--;) {
695 0 : if (((AliAODTrack*)GetTrack(iTrack))->IsMuonGlobalTrack()) {
696 0 : nMuonGlobalTracks++;
697 0 : }
698 : }
699 :
700 0 : return nMuonGlobalTracks;
701 0 : }
702 :
703 : //______________________________________________________________________________
704 : void AliAODEvent::ReadFromTree(TTree *tree, Option_t* opt /*= ""*/)
705 : {
706 : // Connects aod event to tree
707 :
708 0 : if(!tree){
709 0 : AliWarning("Zero Pointer to Tree \n");
710 0 : return;
711 : }
712 : // load the TTree
713 0 : if(!tree->GetTree())tree->LoadTree(0);
714 :
715 : // Try to find AliAODEvent
716 : AliAODEvent *aodEvent = 0;
717 0 : aodEvent = (AliAODEvent*)tree->GetTree()->GetUserInfo()->FindObject("AliAODEvent");
718 0 : if(aodEvent){
719 : // This event is connected to the tree by definition, just say so
720 0 : aodEvent->SetConnected();
721 : // Check if already connected to tree
722 0 : TList* connectedList = (TList*) (tree->GetUserInfo()->FindObject("AODObjectsConnectedToTree"));
723 0 : if (connectedList && (!strcmp(opt, "reconnect"))) {
724 : // If connected use the connected list of objects
725 0 : if (fAODObjects != connectedList) {
726 0 : delete fAODObjects;
727 0 : fAODObjects = connectedList;
728 0 : }
729 0 : GetStdContent();
730 0 : fConnected = kTRUE;
731 0 : return;
732 : }
733 : // Connect to tree
734 : // prevent a memory leak when reading back the TList
735 : // if (!(strcmp(opt, "reconnect"))) fAODObjects->Delete();
736 :
737 : // create a new TList from the UserInfo TList...
738 : // copy constructor does not work...
739 : // fAODObjects = (TList*)(aodEvent->GetList()->Clone());
740 0 : fAODObjects = (TList*)aodEvent->GetList();
741 0 : fAODObjects->SetOwner(kTRUE);
742 0 : if(fAODObjects->GetEntries()<kAODListN)
743 : {
744 0 : AliWarning(Form("AliAODEvent::ReadFromTree() TList contains less than the standard contents %d < %d"
745 : " That might be fine though (at least for filtered AODs)",fAODObjects->GetEntries(),kAODListN));
746 0 : }
747 : //
748 : // Let's find out whether we have friends
749 0 : TList* friendL = tree->GetTree()->GetListOfFriends();
750 0 : if (friendL)
751 : {
752 0 : TIter next(friendL);
753 : TFriendElement* fe;
754 0 : while ((fe = (TFriendElement*)next())){
755 0 : aodEvent = (AliAODEvent*)(fe->GetTree()->GetUserInfo()->FindObject("AliAODEvent"));
756 0 : if (!aodEvent) {
757 0 : printf("No UserInfo on tree \n");
758 : } else {
759 :
760 : // TList* objL = (TList*)(aodEvent->GetList()->Clone());
761 0 : TList* objL = (TList*)aodEvent->GetList();
762 0 : printf("Get list of object from tree %d !!\n", objL->GetEntries());
763 0 : TIter nextobject(objL);
764 : TObject* obj = 0;
765 0 : while((obj = nextobject()))
766 : {
767 0 : printf("Adding object from friend %s !\n", obj->GetName());
768 0 : fAODObjects->Add(obj);
769 : } // object "branch" loop
770 0 : } // has userinfo
771 : } // friend loop
772 0 : } // has friends
773 : // set the branch addresses
774 0 : TIter next(fAODObjects);
775 : TNamed *el;
776 0 : while((el=(TNamed*)next())){
777 0 : TString bname(el->GetName());
778 : // check if branch exists under this Name
779 0 : TBranch *br = tree->GetTree()->GetBranch(bname.Data());
780 0 : if(br){
781 0 : tree->SetBranchAddress(bname.Data(),fAODObjects->GetObjectRef(el));
782 : } else {
783 0 : br = tree->GetBranch(Form("%s.",bname.Data()));
784 0 : if(br){
785 0 : tree->SetBranchAddress(Form("%s.",bname.Data()),fAODObjects->GetObjectRef(el));
786 : }
787 : else{
788 0 : printf("%s %d AliAODEvent::ReadFromTree() No Branch found with Name %s. \n",
789 0 : (char*)__FILE__,__LINE__,bname.Data());
790 : }
791 : }
792 0 : }
793 0 : GetStdContent();
794 : // when reading back we are not owner of the list
795 : // must not delete it
796 0 : fAODObjects->SetOwner(kTRUE);
797 0 : fAODObjects->SetName("AODObjectsConnectedToTree");
798 : // we are not owner of the list objects
799 : // must not delete it
800 0 : tree->GetUserInfo()->Add(fAODObjects);
801 0 : fConnected = kTRUE;
802 0 : }// no aodEvent
803 : else {
804 : // we can't get the list from the user data, create standard content
805 : // and set it by hand
806 0 : CreateStdContent();
807 0 : TIter next(fAODObjects);
808 : TNamed *el;
809 0 : while((el=(TNamed*)next())){
810 0 : TString bname(el->GetName());
811 0 : tree->SetBranchAddress(bname.Data(),fAODObjects->GetObjectRef(el));
812 0 : }
813 0 : GetStdContent();
814 : // when reading back we are not owner of the list
815 : // must not delete it
816 0 : fAODObjects->SetOwner(kTRUE);
817 0 : }
818 0 : }
819 : //______________________________________________________________________________
820 : Int_t AliAODEvent::GetNumberOfPileupVerticesSPD() const{
821 : // count number of SPD pileup vertices
822 0 : Int_t nVertices=GetNumberOfVertices();
823 : Int_t nPileupVertices=0;
824 0 : for(Int_t iVert=0; iVert<nVertices; iVert++){
825 0 : AliAODVertex *v=GetVertex(iVert);
826 0 : if(v->GetType()==AliAODVertex::kPileupSPD) nPileupVertices++;
827 : }
828 0 : return nPileupVertices;
829 : }
830 : //______________________________________________________________________________
831 : Int_t AliAODEvent::GetNumberOfPileupVerticesTracks() const{
832 : // count number of track pileup vertices
833 0 : Int_t nVertices=GetNumberOfVertices();
834 : Int_t nPileupVertices=0;
835 0 : for(Int_t iVert=0; iVert<nVertices; iVert++){
836 0 : AliAODVertex *v=GetVertex(iVert);
837 0 : if(v->GetType()==AliAODVertex::kPileupTracks) nPileupVertices++;
838 : }
839 0 : return nPileupVertices;
840 : }
841 : //______________________________________________________________________________
842 : AliAODVertex* AliAODEvent::GetPrimaryVertexSPD() const{
843 : // Get SPD primary vertex
844 0 : Int_t nVertices=GetNumberOfVertices();
845 0 : for(Int_t iVert=0; iVert<nVertices; iVert++){
846 0 : AliAODVertex *v=GetVertex(iVert);
847 0 : if(v->GetType()==AliAODVertex::kMainSPD) return v;
848 0 : }
849 0 : return 0;
850 0 : }
851 : //______________________________________________________________________________
852 : AliAODVertex* AliAODEvent::GetPrimaryVertexTPC() const{
853 : // Get SPD primary vertex
854 0 : Int_t nVertices=GetNumberOfVertices();
855 0 : for(Int_t iVert=0; iVert<nVertices; iVert++){
856 0 : AliAODVertex *v=GetVertex(iVert);
857 0 : if(v->GetType()==AliAODVertex::kMainTPC) return v;
858 0 : }
859 0 : return 0;
860 0 : }
861 : //______________________________________________________________________________
862 : AliAODVertex* AliAODEvent::GetPileupVertexSPD(Int_t iV) const{
863 : // Get pile-up vertex iV
864 0 : Int_t nVertices=GetNumberOfVertices();
865 : Int_t counter=0;
866 0 : for(Int_t iVert=0; iVert<nVertices; iVert++){
867 0 : AliAODVertex *v=GetVertex(iVert);
868 0 : if(v->GetType()==AliAODVertex::kPileupSPD){
869 0 : if(counter==iV) return v;
870 0 : ++counter;
871 0 : }
872 0 : }
873 0 : return 0;
874 0 : }
875 : //______________________________________________________________________________
876 : AliAODVertex* AliAODEvent::GetPileupVertexTracks(Int_t iV) const{
877 : // Get pile-up vertex iV
878 0 : Int_t nVertices=GetNumberOfVertices();
879 : Int_t counter=0;
880 0 : for(Int_t iVert=0; iVert<nVertices; iVert++){
881 0 : AliAODVertex *v=GetVertex(iVert);
882 0 : if(v->GetType()==AliAODVertex::kPileupTracks){
883 0 : if(counter==iV) return v;
884 0 : ++counter;
885 0 : }
886 0 : }
887 0 : return 0;
888 0 : }
889 : //______________________________________________________________________________
890 : Bool_t AliAODEvent::IsPileupFromSPD(Int_t minContributors,
891 : Double_t minZdist,
892 : Double_t nSigmaZdist,
893 : Double_t nSigmaDiamXY,
894 : Double_t nSigmaDiamZ) const{
895 : //
896 : // This function checks if there was a pile up
897 : // reconstructed with SPD
898 : //
899 0 : AliAODVertex *mainV=GetPrimaryVertexSPD();
900 0 : if(!mainV) return kFALSE;
901 0 : Int_t nc1=mainV->GetNContributors();
902 0 : if(nc1<1) return kFALSE;
903 0 : Int_t nPileVert=GetNumberOfPileupVerticesSPD();
904 0 : if(nPileVert==0) return kFALSE;
905 0 : Int_t nVertices=GetNumberOfVertices();
906 :
907 0 : for(Int_t iVert=0; iVert<nVertices; iVert++){
908 0 : AliAODVertex *pv=GetVertex(iVert);
909 0 : if(pv->GetType()!=AliAODVertex::kPileupSPD) continue;
910 0 : Int_t nc2=pv->GetNContributors();
911 0 : if(nc2>=minContributors){
912 0 : Double_t z1=mainV->GetZ();
913 0 : Double_t z2=pv->GetZ();
914 0 : Double_t distZ=TMath::Abs(z2-z1);
915 0 : Double_t distZdiam=TMath::Abs(z2-GetDiamondZ());
916 0 : Double_t cutZdiam=nSigmaDiamZ*TMath::Sqrt(GetSigma2DiamondZ());
917 0 : if(GetSigma2DiamondZ()<0.0001)cutZdiam=99999.; //protection for missing z diamond information
918 0 : if(distZ>minZdist && distZdiam<cutZdiam){
919 0 : Double_t x2=pv->GetX();
920 0 : Double_t y2=pv->GetY();
921 0 : Double_t distXdiam=TMath::Abs(x2-GetDiamondX());
922 0 : Double_t distYdiam=TMath::Abs(y2-GetDiamondY());
923 0 : Double_t cov1[6],cov2[6];
924 0 : mainV->GetCovarianceMatrix(cov1);
925 0 : pv->GetCovarianceMatrix(cov2);
926 0 : Double_t errxDist=TMath::Sqrt(cov2[0]+GetSigma2DiamondX());
927 0 : Double_t erryDist=TMath::Sqrt(cov2[2]+GetSigma2DiamondY());
928 0 : Double_t errzDist=TMath::Sqrt(cov1[5]+cov2[5]);
929 0 : Double_t cutXdiam=nSigmaDiamXY*errxDist;
930 0 : if(GetSigma2DiamondX()<0.0001)cutXdiam=99999.; //protection for missing diamond information
931 0 : Double_t cutYdiam=nSigmaDiamXY*erryDist;
932 0 : if(GetSigma2DiamondY()<0.0001)cutYdiam=99999.; //protection for missing diamond information
933 0 : if( (distXdiam<cutXdiam) && (distYdiam<cutYdiam) && (distZ>nSigmaZdist*errzDist) ){
934 0 : return kTRUE;
935 : }
936 0 : }
937 0 : }
938 0 : }
939 0 : return kFALSE;
940 0 : }
941 :
942 : //______________________________________________________________________________
943 : void AliAODEvent::Print(Option_t *) const
944 : {
945 : // Print the names of the all branches
946 0 : TIter next(fAODObjects);
947 : TNamed *el;
948 0 : Printf(">>>>> AOD Content <<<<<");
949 0 : while((el=(TNamed*)next())){
950 0 : Printf(">> %s ",el->GetName());
951 : }
952 0 : Printf(">>>>> <<<<<");
953 :
954 : return;
955 0 : }
956 :
957 : //______________________________________________________________________________
958 : void AliAODEvent::AssignIDtoCollection(const TCollection* col)
959 : {
960 : // Static method which assigns a ID to each object in a collection
961 : // In this way the objects are marked as referenced and written with
962 : // an ID. This has the advantage that TRefs to this objects can be
963 : // written by a subsequent process.
964 184 : TIter next(col);
965 : TObject* obj;
966 6206 : while ((obj = next()))
967 2965 : TProcessID::AssignID(obj);
968 92 : }
969 :
970 : //______________________________________________________________________________
971 : Bool_t AliAODEvent::IsPileupFromSPDInMultBins() const {
972 0 : Int_t nTracklets=GetTracklets()->GetNumberOfTracklets();
973 0 : if(nTracklets<20) return IsPileupFromSPD(3,0.8);
974 0 : else if(nTracklets<50) return IsPileupFromSPD(4,0.8);
975 0 : else return IsPileupFromSPD(5,0.8);
976 0 : }
977 :
978 : //______________________________________________________________________________
979 : void AliAODEvent::Reset()
980 : {
981 : // Handle the cases
982 : // Std content + Non std content
983 :
984 0 : ClearStd();
985 :
986 0 : if(fAODObjects->GetSize()>kAODListN){
987 : // we have non std content
988 : // this also covers aodfriends
989 0 : for(int i = kAODListN;i < fAODObjects->GetSize();++i){
990 0 : TObject *pObject = fAODObjects->At(i);
991 : // TClonesArrays
992 0 : if(pObject->InheritsFrom(TClonesArray::Class())){
993 0 : ((TClonesArray*)pObject)->Delete();
994 0 : }
995 0 : else if(!pObject->InheritsFrom(TCollection::Class())){
996 0 : TClass *pClass = TClass::GetClass(pObject->ClassName());
997 0 : if (pClass && pClass->GetListOfMethods()->FindObject("Clear")) {
998 0 : AliDebug(1, Form("Clear for object %s class %s", pObject->GetName(), pObject->ClassName()));
999 0 : pObject->Clear();
1000 0 : }
1001 : else {
1002 0 : AliDebug(1, Form("ResetWithPlacementNew for object %s class %s", pObject->GetName(), pObject->ClassName()));
1003 0 : Long_t dtoronly = TObject::GetDtorOnly();
1004 0 : TObject::SetDtorOnly(pObject);
1005 0 : delete pObject;
1006 0 : pClass->New(pObject);
1007 0 : TObject::SetDtorOnly((void*)dtoronly);
1008 : }
1009 0 : }
1010 : else{
1011 0 : AliWarning(Form("No reset for %s \n",
1012 : pObject->ClassName()));
1013 : }
1014 : }
1015 0 : }
1016 0 : }
1017 :
1018 : // FIXME: Why is this in event and not in header?
1019 : Float_t AliAODEvent::GetVZEROEqMultiplicity(Int_t i) const
1020 : {
1021 : // Get VZERO Multiplicity for channel i
1022 : // Themethod uses the equalization factors
1023 : // stored in the ESD-run object in order to
1024 : // get equal multiplicities within a VZERO rins (1/8 of VZERO)
1025 0 : if (!fAODVZERO || !fHeader) return -1;
1026 :
1027 0 : Int_t ring = i/8;
1028 : Float_t factorSum = 0;
1029 0 : for(Int_t j = 8*ring; j < (8*ring+8); ++j) {
1030 0 : factorSum += fHeader->GetVZEROEqFactors(j);
1031 : }
1032 0 : Float_t factor = fHeader->GetVZEROEqFactors(i)*8./factorSum;
1033 :
1034 0 : return (fAODVZERO->GetMultiplicity(i)/factor);
1035 0 : }
1036 :
1037 : //------------------------------------------------------------
1038 : void AliAODEvent::SetTOFHeader(const AliTOFHeader *header)
1039 : {
1040 : //
1041 : // Set the TOF event_time
1042 : //
1043 :
1044 16 : if (fTOFHeader) {
1045 8 : *fTOFHeader=*header;
1046 : //fTOFHeader->SetName(fgkESDListName[kTOFHeader]);
1047 8 : }
1048 : else {
1049 : // for analysis of reconstructed events
1050 : // when this information is not avaliable
1051 0 : fTOFHeader = new AliTOFHeader(*header);
1052 : //AddObject(fTOFHeader);
1053 : }
1054 :
1055 8 : }
1056 : //------------------------------------------------------------
1057 : AliAODHMPIDrings *AliAODEvent::GetHMPIDringForTrackID(Int_t trackID) const
1058 : {
1059 : //
1060 : // Returns the HMPID object if any for a given track ID
1061 : //
1062 0 : if(GetHMPIDrings())
1063 : {
1064 0 : for(Int_t ien = 0 ; ien < GetNHMPIDrings(); ien++)
1065 : {
1066 0 : if( GetHMPIDring(ien)->GetHmpTrkID() == trackID ) return GetHMPIDring(ien);
1067 : }//rings loop
1068 : }
1069 0 : return 0;
1070 0 : }
1071 : //------------------------------------------------------------
1072 : Int_t AliAODEvent::GetNHMPIDrings() const
1073 : {
1074 : //
1075 : // If there is a list of HMPID rings in the given AOD event, return their number
1076 : //
1077 0 : if ( fHMPIDrings) return fHMPIDrings->GetEntriesFast();
1078 0 : else return -1;
1079 0 : }
1080 : //------------------------------------------------------------
1081 : AliAODHMPIDrings *AliAODEvent::GetHMPIDring(Int_t nRings) const
1082 : {
1083 : //
1084 : // If there is a list of HMPID rings in the given AOD event, return corresponding ring
1085 : //
1086 0 : if(fHMPIDrings) {
1087 0 : if( (AliAODHMPIDrings*)fHMPIDrings->UncheckedAt(nRings) ) {
1088 0 : return (AliAODHMPIDrings*)fHMPIDrings->UncheckedAt(nRings);
1089 : }
1090 0 : else return 0x0;
1091 : }
1092 0 : else return 0x0;
1093 0 : }
1094 : //------------------------------------------------------------
1095 : AliAODTrdTrack& AliAODEvent::AddTrdTrack(const AliVTrdTrack *track) {
1096 6957 : return *(new ((*fTrdTracks)[fTrdTracks->GetEntriesFast()]) AliAODTrdTrack(*track));
1097 0 : }
1098 :
1099 : //______________________________________________________________________________
1100 : void AliAODEvent::ConnectTracks() {
1101 : // Connect tracks to this event
1102 32 : if (fTracksConnected || !fTracks || !fTracks->GetEntriesFast()) return;
1103 : AliAODTrack *track = 0;
1104 24 : track = dynamic_cast<AliAODTrack*>(GetTrack(0));
1105 8 : if(!track) {
1106 0 : AliWarning("Not an AliAODTrack, this is not a standard AOD");
1107 0 : return;
1108 : }
1109 :
1110 8 : TIter next(fTracks);
1111 387 : while ((track=(AliAODTrack*)next())) track->SetAODEvent(this);
1112 8 : fTracksConnected = kTRUE;
1113 16 : }
1114 :
1115 : //______________________________________________________________________________
1116 : Bool_t AliAODEvent::IsIncompleteDAQ()
1117 : {
1118 : // check if DAQ has set the incomplete event attributes
1119 0 : UInt_t daqAttr = GetDAQAttributes();
1120 0 : return (daqAttr&ATTR_2_B(ATTR_INCOMPLETE_EVENT))!=0
1121 0 : || (daqAttr&ATTR_2_B(ATTR_FLUSHED_EVENT))!=0;
1122 :
1123 : }
1124 :
1125 0 : AliVEvent::EDataLayoutType AliAODEvent::GetDataLayoutType() const {return AliVEvent::kAOD;}
1126 :
1127 : //______________________________________________________________________________
1128 : TClonesArray* AliAODEvent::GetDimuons() const
1129 : {
1130 0 : return fDimuons;
1131 : }
1132 :
1133 : //______________________________________________________________________________
1134 : Int_t AliAODEvent::GetNDimuons() const
1135 : {
1136 0 : return fDimuons ? fDimuons->GetEntriesFast() : 0;
1137 : }
1138 :
1139 : //______________________________________________________________________________
1140 : Int_t AliAODEvent::GetNumberOfDimuons() const
1141 : {
1142 0 : return GetNDimuons();
1143 : }
1144 :
1145 : //______________________________________________________________________________
1146 : AliAODDimuon* AliAODEvent::GetDimuon(Int_t nDimu) const
1147 : {
1148 0 : if ( fDimuons )
1149 : {
1150 0 : return (AliAODDimuon*)fDimuons->UncheckedAt(nDimu);
1151 : }
1152 0 : return 0x0;
1153 0 : }
1154 :
1155 : //______________________________________________________________________________
1156 : Int_t AliAODEvent::AddDimuon(const AliAODDimuon* dimu)
1157 : {
1158 0 : if ( fDimuons )
1159 : {
1160 0 : new((*fDimuons)[fDimuons->GetEntriesFast()]) AliAODDimuon(*dimu);
1161 0 : return fDimuons->GetEntriesFast()-1;
1162 : }
1163 0 : return -1;
1164 0 : }
|