]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliAODEvent.cxx
Setter added.
[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 <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
35 ClassImp(AliAODEvent)
36
37 // definition of std AOD member names
38   const char* AliAODEvent::fAODListName[kAODListN] = {"header",
39                                                       "tracks",
40                                                       "vertices",
41                                                       "v0s",
42                                                       "cascades",
43                                                       "tracklets",
44                                                       "jets",
45                                                       "emcalCells",
46                                                       "phosCells",
47                                                       "caloClusters",
48                                                       "fmdClusters",
49                                                       "pmdClusters",
50                                                       "dimuons",
51                                                       "AliAODVZERO"
52                                                       
53 };
54 //______________________________________________________________________________
55 AliAODEvent::AliAODEvent() :
56   AliVEvent(),
57   fAODObjects(new TList()),
58   fAODFolder(0),
59   fConnected(kFALSE),
60   fHeader(0),
61   fTracks(0),
62   fVertices(0),
63   fV0s(0),
64   fCascades(0),
65   fTracklets(0),
66   fJets(0),
67   fEmcalCells(0),
68   fPhosCells(0),
69   fCaloClusters(0),
70   fFmdClusters(0),
71   fPmdClusters(0),
72   fDimuons(0),
73   fAODVZERO(0)
74 {
75   // default constructor
76 }
77
78 //______________________________________________________________________________
79 AliAODEvent::AliAODEvent(const AliAODEvent& aod):
80   AliVEvent(aod),
81   fAODObjects(new TList()),
82   fAODFolder(0),
83   fConnected(kFALSE),
84   fHeader(new AliAODHeader(*aod.fHeader)),
85   fTracks(new TClonesArray(*aod.fTracks)),
86   fVertices(new TClonesArray(*aod.fVertices)),
87   fV0s(new TClonesArray(*aod.fV0s)),
88   fCascades(new TClonesArray(*aod.fCascades)),
89   fTracklets(new AliAODTracklets(*aod.fTracklets)),
90   fJets(new TClonesArray(*aod.fJets)),
91   fEmcalCells(new AliAODCaloCells(*aod.fEmcalCells)),
92   fPhosCells(new AliAODCaloCells(*aod.fPhosCells)),
93   fCaloClusters(new TClonesArray(*aod.fCaloClusters)),
94   fFmdClusters(new TClonesArray(*aod.fFmdClusters)),
95   fPmdClusters(new TClonesArray(*aod.fPmdClusters)),
96   fDimuons(new TClonesArray(*aod.fDimuons)),
97   fAODVZERO(new AliAODVZERO(*aod.fAODVZERO))
98 {
99   // Copy constructor
100   AddObject(fHeader);
101   AddObject(fTracks);
102   AddObject(fVertices);
103   AddObject(fV0s);
104   AddObject(fCascades);
105   AddObject(fTracklets);
106   AddObject(fJets);
107   AddObject(fEmcalCells);
108   AddObject(fPhosCells);
109   AddObject(fCaloClusters);
110   AddObject(fFmdClusters);
111   AddObject(fPmdClusters);
112   AddObject(fDimuons);
113   AddObject(fAODVZERO);
114   fConnected = aod.fConnected;
115   GetStdContent();
116   CreateStdFolders();
117 }
118
119 //______________________________________________________________________________
120 AliAODEvent & AliAODEvent::operator=(const AliAODEvent& aod) {
121
122     // Assignment operator
123
124   if(&aod == this) return *this;
125   AliVEvent::operator=(aod);
126
127   // This assumes that the list is already created
128   // and that the virtual void Copy(Tobject&) function
129   // is correctly implemented in the derived class
130   // otherwise only TObject::Copy() will be used
131
132   if((fAODObjects->GetSize()==0)&&(aod.fAODObjects->GetSize()>=kAODListN)){
133     // We cover the case that we do not yet have the 
134     // standard content but the source has it
135     CreateStdContent();
136   }
137   
138   // Here we have the standard content without user additions, but the content is 
139   // not matching the aod source.
140   
141   // Iterate the list of source objects
142   TIter next(aod.GetList());
143   TObject *its = 0;
144   TString name;
145   while ((its = next())) {
146     name = its->GetName();
147     // Check if we have this object type in out list
148     TObject *mine = fAODObjects->FindObject(name);    
149     if(!mine) {
150       // We have to create the same type of object.
151       TClass* pClass=TClass::GetClass(its->ClassName());     
152       if (!pClass) {
153         AliWarning(Form("Can not find class description for entry %s (%s)\n",
154                    its->ClassName(), name.Data()));
155         continue;
156       }
157       mine=(TObject*)pClass->New();
158       if(!mine){
159         // not in this: can be added to list
160         AliWarning(Form("%s:%d Could not find %s for copying \n",
161                    (char*)__FILE__,__LINE__,name.Data()));
162         continue;
163       }  
164       if(mine->InheritsFrom("TNamed"))  {
165         ((TNamed*)mine)->SetName(name);
166       } else if(mine->InheritsFrom("TCollection")){
167         if(mine->InheritsFrom("TClonesArray")) {
168           TClonesArray *itscl = dynamic_cast<TClonesArray*>(its);
169           if (!itscl) {
170             AliWarning(Form("Class description for entry %s (%s) not TClonesArray\n",
171                    its->ClassName(), name.Data()));
172             continue;
173           
174           }
175                dynamic_cast<TClonesArray*>(mine)->SetClass(itscl->GetClass(), itscl->GetSize());
176         }
177         dynamic_cast<TCollection*>(mine)->SetName(name);
178       }
179       AliDebug(1, Form("adding object %s of type %s", mine->GetName(), mine->ClassName()));
180       AddObject(mine);
181     }
182     // Now we have an object of the same type and name, but different content.        
183     if(!its->InheritsFrom("TCollection")){
184       // simple objects (do they have a Copy method that calls operator= ?)
185       its->Copy(*mine);
186     } else if (its->InheritsFrom("TClonesArray")) {
187       // Create or expand the tclonesarray pointers
188       // so we can directly copy to the object
189       TClonesArray *its_tca = (TClonesArray*)its;
190       TClonesArray *mine_tca = (TClonesArray*)mine;
191       // this leaves the capacity of the TClonesArray the same
192       // except for a factor of 2 increase when size > capacity
193       // does not release any memory occupied by the tca
194       Int_t its_entries = its_tca->GetEntriesFast();
195       mine_tca->ExpandCreate(its_entries);
196       for(int i=0; i<its_entries; i++){
197         // copy 
198         TObject *mine_tca_obj = mine_tca->At(i);
199         TObject *its_tca_obj = its_tca->At(i);
200         // no need to delete first
201         // pointers within the class should be handled by Copy()...
202         // Can there be Empty slots?
203         its_tca_obj->Copy(*mine_tca_obj);
204       }
205     } else {
206       AliWarning(Form("%s:%d cannot copy TCollection \n",
207                       (char*)__FILE__,__LINE__));
208     }
209   }  
210   fConnected = aod.fConnected;
211   return *this;
212 }
213
214
215 //______________________________________________________________________________
216 AliAODEvent::~AliAODEvent() 
217 {
218 // destructor
219     if(fAODObjects&&!fConnected)
220     {
221         delete fAODObjects;
222         fAODObjects = 0;
223     }
224
225     delete fAODFolder;
226     fAODFolder = 0;
227 }
228
229 //______________________________________________________________________________
230 void AliAODEvent::AddObject(TObject* obj) 
231 {
232   // Add an object to the list of objects.
233   // Please be aware that in order to increase performance you should
234   // refrain from using TObjArrays (if possible). Use TClonesArrays, instead.
235   
236   if ( !fAODObjects->FindObject(obj) ) 
237   {
238     fAODObjects->AddLast(obj);
239   }
240 }
241
242 //______________________________________________________________________________
243 void AliAODEvent::RemoveObject(TObject* obj) 
244 {
245   // Removes an object from the list of objects.
246   
247   fAODObjects->Remove(obj);
248 }
249
250 //______________________________________________________________________________
251 TObject *AliAODEvent::FindListObject(const char *objName) const
252 {
253   // Return the pointer to the object with the given name.
254
255   return fAODObjects->FindObject(objName);
256 }
257
258 //______________________________________________________________________________
259 void AliAODEvent::CreateStdContent() 
260 {
261   // create the standard AOD content and set pointers
262
263   // create standard objects and add them to the TList of objects
264   AddObject(new AliAODHeader());
265   AddObject(new TClonesArray("AliAODTrack", 0));
266   AddObject(new TClonesArray("AliAODVertex", 0));
267   AddObject(new TClonesArray("AliAODv0", 0));
268   AddObject(new TClonesArray("AliAODcascade", 0));
269   AddObject(new AliAODTracklets());
270   AddObject(new TClonesArray("AliAODJet", 0));
271   AddObject(new AliAODCaloCells());
272   AddObject(new AliAODCaloCells());
273   AddObject(new TClonesArray("AliAODCaloCluster", 0));
274   AddObject(new TClonesArray("AliAODFmdCluster", 0));
275   AddObject(new TClonesArray("AliAODPmdCluster", 0));
276   AddObject(new TClonesArray("AliAODDimuon", 0));
277   AddObject(new AliAODVZERO());
278   // set names
279   SetStdNames();
280
281   // read back pointers
282   GetStdContent();
283   CreateStdFolders();
284   return;
285 }
286
287 void  AliAODEvent::MakeEntriesReferencable()
288 {
289     // Make all entries referencable in a subsequent process
290     //
291     TIter next(fAODObjects);
292     TObject* obj;
293     while ((obj = next()))
294     {
295         if(obj->InheritsFrom("TCollection"))
296             {
297                 AssignIDtoCollection((TCollection*)obj);
298             }
299     }
300 }
301
302 //______________________________________________________________________________
303 void AliAODEvent::SetStdNames()
304 {
305   // introduce the standard naming
306
307   if(fAODObjects->GetEntries()==kAODListN){
308     for(int i = 0;i < fAODObjects->GetEntries();i++){
309       TObject *fObj = fAODObjects->At(i);
310       if(fObj->InheritsFrom("TNamed")){
311         ((TNamed*)fObj)->SetName(fAODListName[i]);
312       }
313       else if(fObj->InheritsFrom("TClonesArray")){
314         ((TClonesArray*)fObj)->SetName(fAODListName[i]);
315       }
316     }
317   }
318   else{
319     printf("%s:%d SetStdNames() Wrong number of Std Entries \n",(char*)__FILE__,__LINE__);
320   }
321
322
323 void AliAODEvent::CreateStdFolders()
324 {
325     // Create the standard folder structure
326   if(fAODFolder)delete fAODFolder;
327     fAODFolder = gROOT->GetRootFolder()->AddFolder("AOD", "AOD");
328     if(fAODObjects->GetEntries()==kAODListN){
329         for(int i = 0;i < fAODObjects->GetEntries();i++){
330             TObject *fObj = fAODObjects->At(i);
331             if(fObj->InheritsFrom("TClonesArray")){
332                 fAODFolder->AddFolder(fAODListName[i], fAODListName[i], (TCollection*) fObj);
333             } else {
334                 fAODFolder->AddFolder(fAODListName[i], fAODListName[i], 0);
335             }
336         }
337     }
338     else{
339         printf("%s:%d CreateStdFolders() Wrong number of Std Entries \n",(char*)__FILE__,__LINE__);
340     }
341
342
343 //______________________________________________________________________________
344 void AliAODEvent::GetStdContent()
345 {
346   // set pointers for standard content
347
348   fHeader        = (AliAODHeader*)fAODObjects->FindObject("header");
349   fTracks        = (TClonesArray*)fAODObjects->FindObject("tracks");
350   fVertices      = (TClonesArray*)fAODObjects->FindObject("vertices");
351   fV0s           = (TClonesArray*)fAODObjects->FindObject("v0s");
352   fCascades      = (TClonesArray*)fAODObjects->FindObject("cascades");
353   fTracklets     = (AliAODTracklets*)fAODObjects->FindObject("tracklets");
354   fJets          = (TClonesArray*)fAODObjects->FindObject("jets");
355   fEmcalCells    = (AliAODCaloCells*)fAODObjects->FindObject("emcalCells");
356   fPhosCells     = (AliAODCaloCells*)fAODObjects->FindObject("phosCells");
357   fCaloClusters  = (TClonesArray*)fAODObjects->FindObject("caloClusters");
358   fFmdClusters   = (TClonesArray*)fAODObjects->FindObject("fmdClusters");
359   fPmdClusters   = (TClonesArray*)fAODObjects->FindObject("pmdClusters");
360   fDimuons       = (TClonesArray*)fAODObjects->FindObject("dimuons");
361   fAODVZERO      = (AliAODVZERO*)fAODObjects->FindObject("AliAODVZERO");
362 }
363
364 //______________________________________________________________________________
365 void AliAODEvent::ResetStd(Int_t trkArrSize, 
366                            Int_t vtxArrSize, 
367                            Int_t v0ArrSize,
368                            Int_t cascadeArrSize,
369                            Int_t jetSize, 
370                            Int_t caloClusSize, 
371                            Int_t fmdClusSize, 
372                            Int_t pmdClusSize,
373                            Int_t dimuonArrSize
374                            )
375 {
376   // deletes content of standard arrays and resets size 
377   
378   fTracks->Delete();
379   if (trkArrSize > fTracks->GetSize()) 
380     fTracks->Expand(trkArrSize);
381
382   fVertices->Delete();
383   if (vtxArrSize > fVertices->GetSize()) 
384     fVertices->Expand(vtxArrSize);
385
386   fV0s->Delete();
387   if (v0ArrSize > fV0s->GetSize()) 
388     fV0s->Expand(v0ArrSize);
389   
390   fCascades->Delete();
391   if (cascadeArrSize > fCascades->GetSize()) 
392     fCascades->Expand(cascadeArrSize);
393   
394   fJets->Delete();
395   if (jetSize > fJets->GetSize())
396     fJets->Expand(jetSize);
397
398   fCaloClusters->Delete();
399   if (caloClusSize > fCaloClusters->GetSize()) 
400     fCaloClusters->Expand(caloClusSize);
401
402   fFmdClusters->Delete();
403   if (fmdClusSize > fFmdClusters->GetSize()) 
404     fFmdClusters->Expand(fmdClusSize);
405
406   fPmdClusters->Delete();
407   if (pmdClusSize > fPmdClusters->GetSize()) 
408     fPmdClusters->Expand(pmdClusSize);
409     
410   fDimuons->Delete();
411   if (dimuonArrSize > fDimuons->GetSize()) 
412     fDimuons->Expand(dimuonArrSize);
413
414   // Reset the tracklets
415   fTracklets->DeleteContainer();
416   fPhosCells->DeleteContainer();  
417   fEmcalCells->DeleteContainer();
418
419 }
420
421 void AliAODEvent::ClearStd()
422 {
423   // clears the standard arrays
424   fHeader        ->Clear();
425   fTracks        ->Delete();
426   fVertices      ->Delete();
427   fV0s           ->Delete();
428   fCascades      ->Delete();
429   fTracklets     ->DeleteContainer();
430   fJets          ->Delete();
431   fEmcalCells    ->DeleteContainer();
432   fPhosCells     ->DeleteContainer();
433   fCaloClusters  ->Delete();
434   fFmdClusters   ->Clear();
435   fPmdClusters   ->Clear();
436   fDimuons       ->Clear();
437 }
438
439 //_________________________________________________________________
440 Int_t AliAODEvent::GetPHOSClusters(TRefArray *clusters) const
441 {
442   // fills the provided TRefArray with all found phos clusters
443   
444   clusters->Clear();
445   
446   AliAODCaloCluster *cl = 0;
447   Bool_t first = kTRUE;
448   for (Int_t i = 0; i < GetNumberOfCaloClusters() ; i++) {
449     if ( (cl = GetCaloCluster(i)) ) {
450       if (cl->IsPHOS()){
451         if(first) {
452           new (clusters) TRefArray(TProcessID::GetProcessWithUID(cl)); 
453           first=kFALSE;
454         }
455         clusters->Add(cl);
456         //printf("IsPHOS cluster %d, E %2.3f Size: %d \n",i,cl->E(),clusters->GetEntriesFast());
457       }
458     }
459   }
460   return clusters->GetEntriesFast();
461 }
462
463 //_________________________________________________________________
464 Int_t AliAODEvent::GetEMCALClusters(TRefArray *clusters) const
465 {
466   // fills the provided TRefArray with all found emcal clusters
467
468   clusters->Clear();
469   AliAODCaloCluster *cl = 0;
470   Bool_t first = kTRUE;
471   for (Int_t i = 0; i < GetNumberOfCaloClusters(); i++) {
472     if ( (cl = GetCaloCluster(i)) ) {
473       if (cl->IsEMCAL()){
474         if(first) {
475           new (clusters) TRefArray(TProcessID::GetProcessWithUID(cl)); 
476           first=kFALSE;
477         }
478         clusters->Add(cl);
479         //printf("IsEMCal cluster %d, E %2.3f Size: %d \n",i,cl->E(),clusters->GetEntriesFast());
480       }
481     }
482   }
483   return clusters->GetEntriesFast();
484 }
485
486
487 //______________________________________________________________________________
488 Int_t AliAODEvent::GetMuonTracks(TRefArray *muonTracks) const
489 {
490   // fills the provided TRefArray with all found muon tracks
491
492   muonTracks->Clear();
493
494   AliAODTrack *track = 0;
495   for (Int_t iTrack = 0; iTrack < GetNTracks(); iTrack++) {
496     track = GetTrack(iTrack);
497     if (track->IsMuonTrack()) {
498       muonTracks->Add(track);
499     }
500   }
501   
502   return muonTracks->GetEntriesFast();
503 }
504
505
506 //______________________________________________________________________________
507 Int_t AliAODEvent::GetNumberOfMuonTracks() const
508 {
509   // get number of muon tracks
510   Int_t nMuonTracks=0;
511   for (Int_t iTrack = 0; iTrack < GetNTracks(); iTrack++) {
512     if ((GetTrack(iTrack))->IsMuonTrack()) {
513        nMuonTracks++;
514     }
515   }
516   
517   return nMuonTracks;
518 }
519
520 void AliAODEvent::ReadFromTree(TTree *tree, Option_t* opt /*= ""*/)
521 {
522     // Connects aod event to tree
523   
524   if(!tree){
525     Printf("%s %d AliAODEvent::ReadFromTree() Zero Pointer to Tree \n",(char*)__FILE__,__LINE__);
526     return;
527   }
528     // load the TTree
529   if(!tree->GetTree())tree->LoadTree(0);
530   
531     // Try to find AliAODEvent
532   AliAODEvent *aodEvent = 0;
533   aodEvent = (AliAODEvent*)tree->GetTree()->GetUserInfo()->FindObject("AliAODEvent");
534   if(aodEvent){
535       // Check if already connected to tree
536     TList* connectedList = (TList*) (tree->GetUserInfo()->FindObject("AODObjectsConnectedToTree"));
537     if (connectedList && (strcmp(opt, "reconnect"))) {
538         // If connected use the connected list of objects
539         fAODObjects->Delete();
540         fAODObjects = connectedList;
541         GetStdContent(); 
542         fConnected = kTRUE;
543         return;
544     } 
545       // Connect to tree
546       // prevent a memory leak when reading back the TList
547 //      if (!(strcmp(opt, "reconnect"))) fAODObjects->Delete();
548     
549       // create a new TList from the UserInfo TList... 
550       // copy constructor does not work...
551     fAODObjects = (TList*)(aodEvent->GetList()->Clone());
552     fAODObjects->SetOwner(kTRUE);
553     if(fAODObjects->GetEntries()<kAODListN)
554     {
555       AliWarning(Form("AliAODEvent::ReadFromTree() TList contains less than the standard contents %d < %d"
556                       " That might be fine though (at least for filtered AODs)",fAODObjects->GetEntries(),kAODListN));
557     }
558       //
559       // Let's find out whether we have friends
560     TList* friendL = tree->GetTree()->GetListOfFriends();
561     if (friendL) 
562     {
563       TIter next(friendL);
564       TFriendElement* fe;
565       while ((fe = (TFriendElement*)next())){
566         aodEvent = (AliAODEvent*)(fe->GetTree()->GetUserInfo()->FindObject("AliAODEvent"));
567         if (!aodEvent) {
568           printf("No UserInfo on tree \n");
569         } else {
570           
571           TList* objL = (TList*)(aodEvent->GetList()->Clone());
572           printf("Get list of object from tree %d !!\n", objL->GetEntries());
573           TIter nextobject(objL);
574           TObject* obj =  0;
575           while((obj = nextobject()))
576           {
577             printf("Adding object from friend %s !\n", obj->GetName());
578             fAODObjects->Add(obj);
579           } // object "branch" loop
580         } // has userinfo  
581       } // friend loop
582     } // has friends    
583     
584     
585       // set the branch addresses
586     TIter next(fAODObjects);
587     TNamed *el;
588     while((el=(TNamed*)next())){
589       TString bname(el->GetName());
590         // check if branch exists under this Name
591       TBranch *br = tree->GetTree()->GetBranch(bname.Data());
592       if(br){
593         tree->SetBranchAddress(bname.Data(),fAODObjects->GetObjectRef(el));
594       } else {
595         br = tree->GetBranch(Form("%s.",bname.Data()));
596         if(br){
597           tree->SetBranchAddress(Form("%s.",bname.Data()),fAODObjects->GetObjectRef(el));
598         }
599         else{
600           printf("%s %d AliAODEvent::ReadFromTree() No Branch found with Name %s. \n",
601                  (char*)__FILE__,__LINE__,bname.Data());
602         }       
603       }
604     }
605     GetStdContent();
606       // when reading back we are not owner of the list 
607       // must not delete it
608     fAODObjects->SetOwner(kTRUE);
609     fAODObjects->SetName("AODObjectsConnectedToTree");
610       // we are not owner of the list objects 
611       // must not delete it
612     tree->GetUserInfo()->Add(fAODObjects);
613     fConnected = kTRUE;
614   }// no aodEvent
615   else {
616       // we can't get the list from the user data, create standard content
617       // and set it by hand
618     CreateStdContent();
619     TIter next(fAODObjects);
620     TNamed *el;
621     while((el=(TNamed*)next())){
622       TString bname(el->GetName());    
623       tree->SetBranchAddress(bname.Data(),fAODObjects->GetObjectRef(el));
624     }
625     GetStdContent();
626       // when reading back we are not owner of the list 
627       // must not delete it
628     fAODObjects->SetOwner(kTRUE);
629   }
630 }
631   //______________________________________________________________________________
632 Int_t  AliAODEvent::GetNumberOfPileupVerticesSPD() const{
633   // count number of SPD pileup vertices
634   Int_t nVertices=GetNumberOfVertices();
635   Int_t nPileupVertices=0;
636   for(Int_t iVert=0; iVert<nVertices; iVert++){
637     AliAODVertex *v=GetVertex(iVert);
638     if(v->GetType()==AliAODVertex::kPileupSPD) nPileupVertices++;
639   }
640   return nPileupVertices;
641 }
642 //______________________________________________________________________________
643 Int_t  AliAODEvent::GetNumberOfPileupVerticesTracks() const{
644   // count number of track pileup vertices
645   Int_t nVertices=GetNumberOfVertices();
646   Int_t nPileupVertices=0;
647   for(Int_t iVert=0; iVert<nVertices; iVert++){
648     AliAODVertex *v=GetVertex(iVert);
649     if(v->GetType()==AliAODVertex::kPileupTracks) nPileupVertices++;
650   }
651   return nPileupVertices;
652 }
653 //______________________________________________________________________________
654 AliAODVertex* AliAODEvent::GetPrimaryVertexSPD() const{
655   //
656   Int_t nVertices=GetNumberOfVertices();
657   for(Int_t iVert=0; iVert<nVertices; iVert++){
658     AliAODVertex *v=GetVertex(iVert);
659     if(v->GetType()==AliAODVertex::kMainSPD) return v;
660   }
661   return 0;
662 }
663 //______________________________________________________________________________
664 AliAODVertex* AliAODEvent::GetPileupVertexSPD(Int_t iV) const{
665   //
666   Int_t nVertices=GetNumberOfVertices();
667   Int_t counter=0;
668   for(Int_t iVert=0; iVert<nVertices; iVert++){
669     AliAODVertex *v=GetVertex(iVert);
670     if(v->GetType()==AliAODVertex::kPileupSPD){
671       if(counter==iV) return v;
672       ++counter;
673     }
674   }
675   return 0;
676 }
677 //______________________________________________________________________________
678 AliAODVertex* AliAODEvent::GetPileupVertexTracks(Int_t iV) const{
679   //
680   Int_t nVertices=GetNumberOfVertices();
681   Int_t counter=0;
682   for(Int_t iVert=0; iVert<nVertices; iVert++){
683     AliAODVertex *v=GetVertex(iVert);
684     if(v->GetType()==AliAODVertex::kPileupTracks){
685       if(counter==iV) return v;
686       ++counter;
687     }
688   }
689   return 0;
690 }
691 //______________________________________________________________________________
692 Bool_t  AliAODEvent::IsPileupFromSPD(Int_t minContributors, 
693                                      Double_t minZdist, 
694                                      Double_t nSigmaZdist, 
695                                      Double_t nSigmaDiamXY, 
696                                      Double_t nSigmaDiamZ) const{
697   //
698   // This function checks if there was a pile up
699   // reconstructed with SPD
700   //
701   AliAODVertex *mainV=GetPrimaryVertexSPD();
702   if(!mainV) return kFALSE;
703   Int_t nc1=mainV->GetNContributors();
704   if(nc1<1) return kFALSE;
705   Int_t nPileVert=GetNumberOfPileupVerticesSPD();
706   if(nPileVert==0) return kFALSE;
707   Int_t nVertices=GetNumberOfVertices();
708   
709   for(Int_t iVert=0; iVert<nVertices; iVert++){
710     AliAODVertex *pv=GetVertex(iVert);
711     if(pv->GetType()!=AliAODVertex::kPileupSPD) continue;
712     Int_t nc2=pv->GetNContributors();
713     if(nc2>=minContributors){
714       Double_t z1=mainV->GetZ();
715       Double_t z2=pv->GetZ();
716       Double_t distZ=TMath::Abs(z2-z1);
717       Double_t distZdiam=TMath::Abs(z2-GetDiamondZ());
718       Double_t cutZdiam=nSigmaDiamZ*TMath::Sqrt(GetSigma2DiamondZ());
719       if(GetSigma2DiamondZ()<0.0001)cutZdiam=99999.; //protection for missing z diamond information
720       if(distZ>minZdist && distZdiam<cutZdiam){
721         Double_t x2=pv->GetX();
722         Double_t y2=pv->GetY();
723         Double_t distXdiam=TMath::Abs(x2-GetDiamondX());
724         Double_t distYdiam=TMath::Abs(y2-GetDiamondY());
725         Double_t cov1[6],cov2[6];       
726         mainV->GetCovarianceMatrix(cov1);
727         pv->GetCovarianceMatrix(cov2);
728         Double_t errxDist=TMath::Sqrt(cov2[0]+GetSigma2DiamondX());
729         Double_t erryDist=TMath::Sqrt(cov2[2]+GetSigma2DiamondY());
730         Double_t errzDist=TMath::Sqrt(cov1[5]+cov2[5]);
731         Double_t cutXdiam=nSigmaDiamXY*errxDist;
732         if(GetSigma2DiamondX()<0.0001)cutXdiam=99999.; //protection for missing diamond information
733         Double_t cutYdiam=nSigmaDiamXY*erryDist;
734         if(GetSigma2DiamondY()<0.0001)cutYdiam=99999.; //protection for missing diamond information
735         if( (distXdiam<cutXdiam) && (distYdiam<cutYdiam) && (distZ>nSigmaZdist*errzDist) ){
736           return kTRUE;
737         }
738       }
739     }
740   }
741   return kFALSE;
742 }
743
744 //______________________________________________________________________________
745 void AliAODEvent::Print(Option_t *) const
746 {
747   // Print the names of the all branches
748   TIter next(fAODObjects);
749   TNamed *el;
750   Printf(">>>>>  AOD  Content <<<<<");    
751   while((el=(TNamed*)next())){
752     Printf(">> %s ",el->GetName());      
753   }
754   Printf(">>>>>                <<<<<");    
755   
756   return;
757 }
758
759 void AliAODEvent::AssignIDtoCollection(TCollection* col)
760 {
761     // Static method which assigns a ID to each object in a collection
762     // In this way the objects are marked as referenced and written with 
763     // an ID. This has the advantage that TRefs to this objects can be 
764     // written by a subsequent process.
765     TIter next(col);
766     TObject* obj;
767     while ((obj = next()))
768         TProcessID::AssignID(obj);
769 }
770
771 Bool_t AliAODEvent::IsPileupFromSPDInMultBins() const {
772     Int_t nTracklets=GetTracklets()->GetNumberOfTracklets();
773     if(nTracklets<20) return IsPileupFromSPD(3,0.8);
774     else if(nTracklets<50) return IsPileupFromSPD(4,0.8);
775     else return IsPileupFromSPD(5,0.8);
776 }
777