]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliAODEvent.cxx
Third step towards VZERO AOD: inclusing in the AliAODEvent and in esd-filtering analy...
[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   fAODObjects->AddLast(obj);
237 }
238
239 //______________________________________________________________________________
240 void AliAODEvent::RemoveObject(TObject* obj) 
241 {
242   // Removes an object from the list of objects.
243   
244   fAODObjects->Remove(obj);
245 }
246
247 //______________________________________________________________________________
248 TObject *AliAODEvent::FindListObject(const char *objName) const
249 {
250   // Return the pointer to the object with the given name.
251
252   return fAODObjects->FindObject(objName);
253 }
254
255 //______________________________________________________________________________
256 void AliAODEvent::CreateStdContent() 
257 {
258   // create the standard AOD content and set pointers
259
260   // create standard objects and add them to the TList of objects
261   AddObject(new AliAODHeader());
262   AddObject(new TClonesArray("AliAODTrack", 0));
263   AddObject(new TClonesArray("AliAODVertex", 0));
264   AddObject(new TClonesArray("AliAODv0", 0));
265   AddObject(new TClonesArray("AliAODcascade", 0));
266   AddObject(new AliAODTracklets());
267   AddObject(new TClonesArray("AliAODJet", 0));
268   AddObject(new AliAODCaloCells());
269   AddObject(new AliAODCaloCells());
270   AddObject(new TClonesArray("AliAODCaloCluster", 0));
271   AddObject(new TClonesArray("AliAODFmdCluster", 0));
272   AddObject(new TClonesArray("AliAODPmdCluster", 0));
273   AddObject(new TClonesArray("AliAODDimuon", 0));
274   AddObject(new AliAODVZERO());
275   // set names
276   SetStdNames();
277
278   // read back pointers
279   GetStdContent();
280   CreateStdFolders();
281   return;
282 }
283
284 void  AliAODEvent::MakeEntriesReferencable()
285 {
286     // Make all entries referencable in a subsequent process
287     //
288     TIter next(fAODObjects);
289     TObject* obj;
290     while ((obj = next()))
291     {
292         if(obj->InheritsFrom("TCollection"))
293             {
294                 AssignIDtoCollection((TCollection*)obj);
295             }
296     }
297 }
298
299 //______________________________________________________________________________
300 void AliAODEvent::SetStdNames()
301 {
302   // introduce the standard naming
303
304   if(fAODObjects->GetEntries()==kAODListN){
305     for(int i = 0;i < fAODObjects->GetEntries();i++){
306       TObject *fObj = fAODObjects->At(i);
307       if(fObj->InheritsFrom("TNamed")){
308         ((TNamed*)fObj)->SetName(fAODListName[i]);
309       }
310       else if(fObj->InheritsFrom("TClonesArray")){
311         ((TClonesArray*)fObj)->SetName(fAODListName[i]);
312       }
313     }
314   }
315   else{
316     printf("%s:%d SetStdNames() Wrong number of Std Entries \n",(char*)__FILE__,__LINE__);
317   }
318
319
320 void AliAODEvent::CreateStdFolders()
321 {
322     // Create the standard folder structure
323   if(fAODFolder)delete fAODFolder;
324     fAODFolder = gROOT->GetRootFolder()->AddFolder("AOD", "AOD");
325     if(fAODObjects->GetEntries()==kAODListN){
326         for(int i = 0;i < fAODObjects->GetEntries();i++){
327             TObject *fObj = fAODObjects->At(i);
328             if(fObj->InheritsFrom("TClonesArray")){
329                 fAODFolder->AddFolder(fAODListName[i], fAODListName[i], (TCollection*) fObj);
330             } else {
331                 fAODFolder->AddFolder(fAODListName[i], fAODListName[i], 0);
332             }
333         }
334     }
335     else{
336         printf("%s:%d CreateStdFolders() Wrong number of Std Entries \n",(char*)__FILE__,__LINE__);
337     }
338
339
340 //______________________________________________________________________________
341 void AliAODEvent::GetStdContent()
342 {
343   // set pointers for standard content
344
345   fHeader        = (AliAODHeader*)fAODObjects->FindObject("header");
346   fTracks        = (TClonesArray*)fAODObjects->FindObject("tracks");
347   fVertices      = (TClonesArray*)fAODObjects->FindObject("vertices");
348   fV0s           = (TClonesArray*)fAODObjects->FindObject("v0s");
349   fCascades      = (TClonesArray*)fAODObjects->FindObject("cascades");
350   fTracklets     = (AliAODTracklets*)fAODObjects->FindObject("tracklets");
351   fJets          = (TClonesArray*)fAODObjects->FindObject("jets");
352   fEmcalCells    = (AliAODCaloCells*)fAODObjects->FindObject("emcalCells");
353   fPhosCells     = (AliAODCaloCells*)fAODObjects->FindObject("phosCells");
354   fCaloClusters  = (TClonesArray*)fAODObjects->FindObject("caloClusters");
355   fFmdClusters   = (TClonesArray*)fAODObjects->FindObject("fmdClusters");
356   fPmdClusters   = (TClonesArray*)fAODObjects->FindObject("pmdClusters");
357   fDimuons       = (TClonesArray*)fAODObjects->FindObject("dimuons");
358   fAODVZERO      = (AliAODVZERO*)fAODObjects->FindObject("AliAODVZERO");
359 }
360
361 //______________________________________________________________________________
362 void AliAODEvent::ResetStd(Int_t trkArrSize, 
363                            Int_t vtxArrSize, 
364                            Int_t v0ArrSize,
365                            Int_t cascadeArrSize,
366                            Int_t jetSize, 
367                            Int_t caloClusSize, 
368                            Int_t fmdClusSize, 
369                            Int_t pmdClusSize,
370                            Int_t dimuonArrSize
371                            )
372 {
373   // deletes content of standard arrays and resets size 
374   
375   fTracks->Delete();
376   if (trkArrSize > fTracks->GetSize()) 
377     fTracks->Expand(trkArrSize);
378
379   fVertices->Delete();
380   if (vtxArrSize > fVertices->GetSize()) 
381     fVertices->Expand(vtxArrSize);
382
383   fV0s->Delete();
384   if (v0ArrSize > fV0s->GetSize()) 
385     fV0s->Expand(v0ArrSize);
386   
387   fCascades->Delete();
388   if (cascadeArrSize > fCascades->GetSize()) 
389     fCascades->Expand(cascadeArrSize);
390   
391   fJets->Delete();
392   if (jetSize > fJets->GetSize())
393     fJets->Expand(jetSize);
394
395   fCaloClusters->Delete();
396   if (caloClusSize > fCaloClusters->GetSize()) 
397     fCaloClusters->Expand(caloClusSize);
398
399   fFmdClusters->Delete();
400   if (fmdClusSize > fFmdClusters->GetSize()) 
401     fFmdClusters->Expand(fmdClusSize);
402
403   fPmdClusters->Delete();
404   if (pmdClusSize > fPmdClusters->GetSize()) 
405     fPmdClusters->Expand(pmdClusSize);
406     
407   fDimuons->Delete();
408   if (dimuonArrSize > fDimuons->GetSize()) 
409     fDimuons->Expand(dimuonArrSize);
410
411   // Reset the tracklets
412   fTracklets->DeleteContainer();
413   fPhosCells->DeleteContainer();  
414   fEmcalCells->DeleteContainer();
415
416 }
417
418 void AliAODEvent::ClearStd()
419 {
420   // clears the standard arrays
421   fHeader        ->Clear();
422   fTracks        ->Delete();
423   fVertices      ->Delete();
424   fV0s           ->Delete();
425   fCascades      ->Delete();
426   fTracklets     ->DeleteContainer();
427   fJets          ->Delete();
428   fEmcalCells    ->DeleteContainer();
429   fPhosCells     ->DeleteContainer();
430   fCaloClusters  ->Delete();
431   fFmdClusters   ->Clear();
432   fPmdClusters   ->Clear();
433   fDimuons       ->Clear();
434 }
435
436 //_________________________________________________________________
437 Int_t AliAODEvent::GetPHOSClusters(TRefArray *clusters) const
438 {
439   // fills the provided TRefArray with all found phos clusters
440   
441   clusters->Clear();
442   
443   AliAODCaloCluster *cl = 0;
444   Bool_t first = kTRUE;
445   for (Int_t i = 0; i < GetNumberOfCaloClusters() ; i++) {
446     if ( (cl = GetCaloCluster(i)) ) {
447       if (cl->IsPHOS()){
448         if(first) {
449           new (clusters) TRefArray(TProcessID::GetProcessWithUID(cl)); 
450           first=kFALSE;
451         }
452         clusters->Add(cl);
453         //printf("IsPHOS cluster %d, E %2.3f Size: %d \n",i,cl->E(),clusters->GetEntriesFast());
454       }
455     }
456   }
457   return clusters->GetEntriesFast();
458 }
459
460 //_________________________________________________________________
461 Int_t AliAODEvent::GetEMCALClusters(TRefArray *clusters) const
462 {
463   // fills the provided TRefArray with all found emcal clusters
464
465   clusters->Clear();
466   AliAODCaloCluster *cl = 0;
467   Bool_t first = kTRUE;
468   for (Int_t i = 0; i < GetNumberOfCaloClusters(); i++) {
469     if ( (cl = GetCaloCluster(i)) ) {
470       if (cl->IsEMCAL()){
471         if(first) {
472           new (clusters) TRefArray(TProcessID::GetProcessWithUID(cl)); 
473           first=kFALSE;
474         }
475         clusters->Add(cl);
476         //printf("IsEMCal cluster %d, E %2.3f Size: %d \n",i,cl->E(),clusters->GetEntriesFast());
477       }
478     }
479   }
480   return clusters->GetEntriesFast();
481 }
482
483
484 //______________________________________________________________________________
485 Int_t AliAODEvent::GetMuonTracks(TRefArray *muonTracks) const
486 {
487   // fills the provided TRefArray with all found muon tracks
488
489   muonTracks->Clear();
490
491   AliAODTrack *track = 0;
492   for (Int_t iTrack = 0; iTrack < GetNTracks(); iTrack++) {
493     track = GetTrack(iTrack);
494     if (track->IsMuonTrack()) {
495       muonTracks->Add(track);
496     }
497   }
498   
499   return muonTracks->GetEntriesFast();
500 }
501
502
503 //______________________________________________________________________________
504 Int_t AliAODEvent::GetNumberOfMuonTracks() const
505 {
506   // get number of muon tracks
507   Int_t nMuonTracks=0;
508   for (Int_t iTrack = 0; iTrack < GetNTracks(); iTrack++) {
509     if ((GetTrack(iTrack))->IsMuonTrack()) {
510        nMuonTracks++;
511     }
512   }
513   
514   return nMuonTracks;
515 }
516
517 void AliAODEvent::ReadFromTree(TTree *tree, Option_t* opt /*= ""*/)
518 {
519     // Connects aod event to tree
520   
521   if(!tree){
522     Printf("%s %d AliAODEvent::ReadFromTree() Zero Pointer to Tree \n",(char*)__FILE__,__LINE__);
523     return;
524   }
525     // load the TTree
526   if(!tree->GetTree())tree->LoadTree(0);
527   
528     // Try to find AliAODEvent
529   AliAODEvent *aodEvent = 0;
530   aodEvent = (AliAODEvent*)tree->GetTree()->GetUserInfo()->FindObject("AliAODEvent");
531   if(aodEvent){
532       // Check if already connected to tree
533     TList* connectedList = (TList*) (tree->GetUserInfo()->FindObject("AODObjectsConnectedToTree"));
534     if (connectedList && (strcmp(opt, "reconnect"))) {
535         // If connected use the connected list of objects
536         fAODObjects->Delete();
537         fAODObjects = connectedList;
538         GetStdContent(); 
539         fConnected = kTRUE;
540         return;
541     } 
542       // Connect to tree
543       // prevent a memory leak when reading back the TList
544 //      if (!(strcmp(opt, "reconnect"))) fAODObjects->Delete();
545     
546       // create a new TList from the UserInfo TList... 
547       // copy constructor does not work...
548     fAODObjects = (TList*)(aodEvent->GetList()->Clone());
549     fAODObjects->SetOwner(kTRUE);
550     if(fAODObjects->GetEntries()<kAODListN){
551       printf("%s %d AliAODEvent::ReadFromTree() TList contains less than the standard contents %d < %d \n",
552              (char*)__FILE__,__LINE__,fAODObjects->GetEntries(),kAODListN);
553     }
554       //
555       // Let's find out whether we have friends
556     TList* friendL = tree->GetTree()->GetListOfFriends();
557     if (friendL) 
558     {
559       TIter next(friendL);
560       TFriendElement* fe;
561       while ((fe = (TFriendElement*)next())){
562         aodEvent = (AliAODEvent*)(fe->GetTree()->GetUserInfo()->FindObject("AliAODEvent"));
563         if (!aodEvent) {
564           printf("No UserInfo on tree \n");
565         } else {
566           
567           TList* objL = (TList*)(aodEvent->GetList()->Clone());
568           printf("Get list of object from tree %d !!\n", objL->GetEntries());
569           TIter nextobject(objL);
570           TObject* obj =  0;
571           while((obj = nextobject()))
572           {
573             printf("Adding object from friend %s !\n", obj->GetName());
574             fAODObjects->Add(obj);
575           } // object "branch" loop
576         } // has userinfo  
577       } // friend loop
578     } // has friends    
579     
580     
581       // set the branch addresses
582     TIter next(fAODObjects);
583     TNamed *el;
584     while((el=(TNamed*)next())){
585       TString bname(el->GetName());
586         // check if branch exists under this Name
587       TBranch *br = tree->GetTree()->GetBranch(bname.Data());
588       if(br){
589         tree->SetBranchAddress(bname.Data(),fAODObjects->GetObjectRef(el));
590       } else {
591         br = tree->GetBranch(Form("%s.",bname.Data()));
592         if(br){
593           tree->SetBranchAddress(Form("%s.",bname.Data()),fAODObjects->GetObjectRef(el));
594         }
595         else{
596           printf("%s %d AliAODEvent::ReadFromTree() No Branch found with Name %s. \n",
597                  (char*)__FILE__,__LINE__,bname.Data());
598         }       
599       }
600     }
601     GetStdContent();
602       // when reading back we are not owner of the list 
603       // must not delete it
604     fAODObjects->SetOwner(kTRUE);
605     fAODObjects->SetName("AODObjectsConnectedToTree");
606       // we are not owner of the list objects 
607       // must not delete it
608     tree->GetUserInfo()->Add(fAODObjects);
609     fConnected = kTRUE;
610   }// no aodEvent
611   else {
612       // we can't get the list from the user data, create standard content
613       // and set it by hand
614     CreateStdContent();
615     TIter next(fAODObjects);
616     TNamed *el;
617     while((el=(TNamed*)next())){
618       TString bname(el->GetName());    
619       tree->SetBranchAddress(bname.Data(),fAODObjects->GetObjectRef(el));
620     }
621     GetStdContent();
622       // when reading back we are not owner of the list 
623       // must not delete it
624     fAODObjects->SetOwner(kTRUE);
625   }
626 }
627   //______________________________________________________________________________
628 Int_t  AliAODEvent::GetNumberOfPileupVerticesSPD() const{
629   // count number of SPD pileup vertices
630   Int_t nVertices=GetNumberOfVertices();
631   Int_t nPileupVertices=0;
632   for(Int_t iVert=0; iVert<nVertices; iVert++){
633     AliAODVertex *v=GetVertex(iVert);
634     if(v->GetType()==AliAODVertex::kPileupSPD) nPileupVertices++;
635   }
636   return nPileupVertices;
637 }
638 //______________________________________________________________________________
639 Int_t  AliAODEvent::GetNumberOfPileupVerticesTracks() const{
640   // count number of track pileup vertices
641   Int_t nVertices=GetNumberOfVertices();
642   Int_t nPileupVertices=0;
643   for(Int_t iVert=0; iVert<nVertices; iVert++){
644     AliAODVertex *v=GetVertex(iVert);
645     if(v->GetType()==AliAODVertex::kPileupTracks) nPileupVertices++;
646   }
647   return nPileupVertices;
648 }
649 //______________________________________________________________________________
650 AliAODVertex* AliAODEvent::GetPrimaryVertexSPD() const{
651   //
652   Int_t nVertices=GetNumberOfVertices();
653   for(Int_t iVert=0; iVert<nVertices; iVert++){
654     AliAODVertex *v=GetVertex(iVert);
655     if(v->GetType()==AliAODVertex::kMainSPD) return v;
656   }
657   return 0;
658 }
659 //______________________________________________________________________________
660 AliAODVertex* AliAODEvent::GetPileupVertexSPD(Int_t iV) const{
661   //
662   Int_t nVertices=GetNumberOfVertices();
663   Int_t counter=0;
664   for(Int_t iVert=0; iVert<nVertices; iVert++){
665     AliAODVertex *v=GetVertex(iVert);
666     if(v->GetType()==AliAODVertex::kPileupSPD){
667       if(counter==iV) return v;
668       ++counter;
669     }
670   }
671   return 0;
672 }
673 //______________________________________________________________________________
674 AliAODVertex* AliAODEvent::GetPileupVertexTracks(Int_t iV) const{
675   //
676   Int_t nVertices=GetNumberOfVertices();
677   Int_t counter=0;
678   for(Int_t iVert=0; iVert<nVertices; iVert++){
679     AliAODVertex *v=GetVertex(iVert);
680     if(v->GetType()==AliAODVertex::kPileupTracks){
681       if(counter==iV) return v;
682       ++counter;
683     }
684   }
685   return 0;
686 }
687 //______________________________________________________________________________
688 Bool_t  AliAODEvent::IsPileupFromSPD(Int_t minContributors, 
689                                      Double_t minZdist, 
690                                      Double_t nSigmaZdist, 
691                                      Double_t nSigmaDiamXY, 
692                                      Double_t nSigmaDiamZ) const{
693   //
694   // This function checks if there was a pile up
695   // reconstructed with SPD
696   //
697   AliAODVertex *mainV=GetPrimaryVertexSPD();
698   if(!mainV) return kFALSE;
699   Int_t nc1=mainV->GetNContributors();
700   if(nc1<1) return kFALSE;
701   Int_t nPileVert=GetNumberOfPileupVerticesSPD();
702   if(nPileVert==0) return kFALSE;
703   Int_t nVertices=GetNumberOfVertices();
704   
705   for(Int_t iVert=0; iVert<nVertices; iVert++){
706     AliAODVertex *pv=GetVertex(iVert);
707     if(pv->GetType()!=AliAODVertex::kPileupSPD) continue;
708     Int_t nc2=pv->GetNContributors();
709     if(nc2>=minContributors){
710       Double_t z1=mainV->GetZ();
711       Double_t z2=pv->GetZ();
712       Double_t distZ=TMath::Abs(z2-z1);
713       Double_t distZdiam=TMath::Abs(z2-GetDiamondZ());
714       Double_t cutZdiam=nSigmaDiamZ*TMath::Sqrt(GetSigma2DiamondZ());
715       if(GetSigma2DiamondZ()<0.0001)cutZdiam=99999.; //protection for missing z diamond information
716       if(distZ>minZdist && distZdiam<cutZdiam){
717         Double_t x2=pv->GetX();
718         Double_t y2=pv->GetY();
719         Double_t distXdiam=TMath::Abs(x2-GetDiamondX());
720         Double_t distYdiam=TMath::Abs(y2-GetDiamondY());
721         Double_t cov1[6],cov2[6];       
722         mainV->GetCovarianceMatrix(cov1);
723         pv->GetCovarianceMatrix(cov2);
724         Double_t errxDist=TMath::Sqrt(cov2[0]+GetSigma2DiamondX());
725         Double_t erryDist=TMath::Sqrt(cov2[2]+GetSigma2DiamondY());
726         Double_t errzDist=TMath::Sqrt(cov1[5]+cov2[5]);
727         Double_t cutXdiam=nSigmaDiamXY*errxDist;
728         if(GetSigma2DiamondX()<0.0001)cutXdiam=99999.; //protection for missing diamond information
729         Double_t cutYdiam=nSigmaDiamXY*erryDist;
730         if(GetSigma2DiamondY()<0.0001)cutYdiam=99999.; //protection for missing diamond information
731         if( (distXdiam<cutXdiam) && (distYdiam<cutYdiam) && (distZ>nSigmaZdist*errzDist) ){
732           return kTRUE;
733         }
734       }
735     }
736   }
737   return kFALSE;
738 }
739
740 //______________________________________________________________________________
741 void AliAODEvent::Print(Option_t *) const
742 {
743   // Print the names of the all branches
744   TIter next(fAODObjects);
745   TNamed *el;
746   Printf(">>>>>  AOD  Content <<<<<");    
747   while((el=(TNamed*)next())){
748     Printf(">> %s ",el->GetName());      
749   }
750   Printf(">>>>>                <<<<<");    
751   
752   return;
753 }
754
755 void AliAODEvent::AssignIDtoCollection(TCollection* col)
756 {
757     // Static method which assigns a ID to each object in a collection
758     // In this way the objects are marked as referenced and written with 
759     // an ID. This has the advantage that TRefs to this objects can be 
760     // written by a subsequent process.
761     TIter next(col);
762     TObject* obj;
763     while ((obj = next()))
764         TProcessID::AssignID(obj);
765 }
766
767 Bool_t AliAODEvent::IsPileupFromSPDInMultBins() const {
768     Int_t nTracklets=GetTracklets()->GetNumberOfTracklets();
769     if(nTracklets<20) return IsPileupFromSPD(3,0.8);
770     else if(nTracklets<50) return IsPileupFromSPD(4,0.8);
771     else return IsPileupFromSPD(5,0.8);
772 }
773