]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliAODHandler.cxx
fixes from Laurent for the MC branch in the AOD filters
[u/mrichter/AliRoot.git] / STEER / AliAODHandler.cxx
1
2 /**************************************************************************
3  * Copyright(c) 1998-2007, ALICE Experiment at CERN, All rights reserved. *
4  *                                                                        *
5  * Author: The ALICE Off-line Project.                                    *
6  * Contributors are mentioned in the code where appropriate.              *
7  *                                                                        *
8  * Permission to use, copy, modify and distribute this software and its   *
9  * documentation strictly for non-commercial purposes is hereby granted   *
10  * without fee, provided that the above copyright notice appears in all   *
11  * copies and that both the copyright notice and this permission notice   *
12  * appear in the supporting documentation. The authors make no claims     *
13  * about the suitability of this software for any purpose. It is          *
14  * provided "as is" without express or implied warranty.                  *
15  **************************************************************************/
16
17 /* $Id$ */
18
19 //-------------------------------------------------------------------------
20 //     Implementation of the Virtual Event Handler Interface for AOD
21 //     Author: Andreas Morsch, CERN
22 //-------------------------------------------------------------------------
23
24
25 #include <TTree.h>
26 #include <TFile.h>
27 #include <TString.h>
28 #include <TList.h>
29 #include <TROOT.h>
30
31 #include "AliLog.h"
32 #include "AliAODHandler.h"
33 #include "AliAODEvent.h"
34 #include "AliAODTracklets.h"
35 #include "AliStack.h"
36 #include "AliAODMCParticle.h"
37 #include "AliAODMCHeader.h"
38 #include "AliMCEventHandler.h"
39 #include "AliMCEvent.h"
40 #include "AliGenEventHeader.h"
41 #include "AliGenHijingEventHeader.h"
42 #include "AliGenDPMjetEventHeader.h"
43 #include "AliGenPythiaEventHeader.h"
44 #include "AliGenCocktailEventHeader.h"
45 #include "AliCodeTimer.h"
46 #include "AliAODBranchReplicator.h"
47 #include "Riostream.h"
48
49 ClassImp(AliAODHandler)
50
51 //______________________________________________________________________________
52 AliAODHandler::AliAODHandler() :
53     AliVEventHandler(),
54     fIsStandard(kTRUE),
55     fFillAOD(kTRUE),
56     fFillAODRun(kTRUE),
57     fNeedsHeaderReplication(kFALSE),
58     fNeedsTracksBranchReplication(kFALSE),
59     fNeedsVerticesBranchReplication(kFALSE),
60     fNeedsV0sBranchReplication(kFALSE),
61     fNeedsCascadesBranchReplication(kFALSE),
62     fNeedsTrackletsBranchReplication(kFALSE),
63     fNeedsPMDClustersBranchReplication(kFALSE),
64     fNeedsJetsBranchReplication(kFALSE),
65     fNeedsFMDClustersBranchReplication(kFALSE),
66     fNeedsCaloClustersBranchReplication(kFALSE),
67     fNeedsMCParticlesBranchReplication(kFALSE),
68     fNeedsDimuonsBranchReplication(kFALSE),
69     fAODIsReplicated(kFALSE),
70     fAODEvent(NULL),
71     fMCEventH(NULL),
72     fTreeA(NULL),
73     fFileA(NULL),
74     fFileName(""),
75     fExtensions(NULL),
76     fFilters(NULL)
77 {
78   // default constructor
79 }
80
81 //______________________________________________________________________________
82 AliAODHandler::AliAODHandler(const char* name, const char* title):
83     AliVEventHandler(name, title),
84     fIsStandard(kTRUE),
85     fFillAOD(kTRUE),
86     fFillAODRun(kTRUE),
87     fNeedsHeaderReplication(kFALSE),
88     fNeedsTracksBranchReplication(kFALSE),
89     fNeedsVerticesBranchReplication(kFALSE),
90     fNeedsV0sBranchReplication(kFALSE),
91     fNeedsCascadesBranchReplication(kFALSE),
92     fNeedsTrackletsBranchReplication(kFALSE),
93     fNeedsPMDClustersBranchReplication(kFALSE),
94     fNeedsJetsBranchReplication(kFALSE),
95     fNeedsFMDClustersBranchReplication(kFALSE),
96     fNeedsCaloClustersBranchReplication(kFALSE),
97     fNeedsMCParticlesBranchReplication(kFALSE),
98     fNeedsDimuonsBranchReplication(kFALSE),
99     fAODIsReplicated(kFALSE),
100     fAODEvent(NULL),
101     fMCEventH(NULL),
102     fTreeA(NULL),
103     fFileA(NULL),
104     fFileName(""),
105     fExtensions(NULL),
106     fFilters(NULL)
107 {
108 // Normal constructor.
109 }
110
111 //______________________________________________________________________________
112 AliAODHandler::~AliAODHandler() 
113 {
114  // Destructor.
115   
116   delete fAODEvent;
117
118   if (fFileA)
119   {
120     // is already handled in TerminateIO
121     fFileA->Close();
122     delete fFileA;
123     fTreeA = 0;
124   }
125   delete fTreeA;
126   delete fExtensions;
127   delete fFilters;
128 }
129
130 //______________________________________________________________________________
131 Bool_t AliAODHandler::Init(Option_t* opt)
132 {
133   // Initialize IO
134   //
135   // Create the AODevent object
136     
137   Bool_t createStdAOD = fIsStandard || fFillAOD;
138   if(!fAODEvent && createStdAOD){
139     fAODEvent = new AliAODEvent();
140     if (fIsStandard) 
141       fAODEvent->CreateStdContent();
142   }
143   //
144   // File opening according to execution mode
145   TString option(opt);
146   option.ToLower();
147   if (createStdAOD) {
148     TDirectory *owd = gDirectory;
149     if (option.Contains("proof")) {
150       // proof
151       // Merging via files. Need to access analysis manager via interpreter.
152       gROOT->ProcessLine(Form("AliAnalysisDataContainer *c_common_out = AliAnalysisManager::GetAnalysisManager()->GetCommonOutputContainer();"));
153       gROOT->ProcessLine(Form("AliAnalysisManager::GetAnalysisManager()->OpenProofFile(c_common_out, \"RECREATE\");"));
154       fFileA = gFile;
155     } else {
156       // local and grid
157       fFileA = new TFile(fFileName.Data(), "RECREATE");
158     }
159     CreateTree(1);
160     owd->cd();
161   }  
162   if (fExtensions) {
163      TIter next(fExtensions);
164      AliAODExtension *ext;
165      while ((ext=(AliAODExtension*)next())) ext->Init(option);
166   }   
167   if (fFilters) {
168      TIter nextf(fFilters);
169      AliAODExtension *filteredAOD;
170      while ((filteredAOD=(AliAODExtension*)nextf())) {
171         filteredAOD->SetEvent(fAODEvent);
172         filteredAOD->Init(option);
173      }
174   }   
175   
176   return kTRUE;
177 }
178
179 //______________________________________________________________________________
180 void AliAODHandler::Print(Option_t* opt) const
181 {
182   // Print info about this object
183
184   cout << opt << Form("IsStandard %d filename=%s",fIsStandard,fFileName.Data()) << endl;
185   
186   if ( fExtensions ) 
187   {
188     cout << opt << fExtensions->GetEntries() << " extensions :" << endl;
189     PrintExtensions(*fExtensions);    
190   }
191   if ( fFilters ) 
192   {
193     cout << opt << fFilters->GetEntries() << " filters :" << endl;
194     PrintExtensions(*fFilters);      
195   }
196 }
197
198 //______________________________________________________________________________
199 void AliAODHandler::PrintExtensions(const TObjArray& array) const
200 {
201   // Show the list of aod extensions
202   TIter next(&array);
203   AliAODExtension* ext(0x0);
204   while ( ( ext = static_cast<AliAODExtension*>(next()) ) )
205   {
206     ext->Print("   ");
207   }
208 }
209
210 //______________________________________________________________________________
211 void AliAODHandler::StoreMCParticles(){
212
213   // 
214   // Remap the labels from ESD stack and store
215   // the AODMCParticles, makes only sense if we have
216   // the mcparticles branch
217   // has to be done here since we cannot know in advance 
218   // which particles are needed (e.g. by the tracks etc.)
219   //
220   // Particles have been selected by AliMCEventhanlder->SelectParticle()
221   // To use the MCEventhandler here we need to set it from the outside
222   // can vanish when Handler go to the ANALYSISalice library
223   //
224   // The Branch booking for mcParticles and mcHeader has to happen 
225   // in an external task for now since the AODHandler does not have access
226   // the AnalysisManager. For the same reason the pointer t o the MCEventH
227   // has to passed to the AOD Handler by this task 
228   // (doing this in the steering macro would not work on PROOF)
229
230   if (!fAODEvent) return;
231   TClonesArray *mcarray = (TClonesArray*)fAODEvent->FindListObject(AliAODMCParticle::StdBranchName()); 
232   if(!mcarray)return;
233
234   AliAODMCHeader *mcHeader = (AliAODMCHeader*)fAODEvent->FindListObject(AliAODMCHeader::StdBranchName()); 
235   if(!mcHeader)return;
236
237   // Get the MC Infos.. Handler needs to be set before 
238   // while adding the branch
239   // This needs to be done, not to depend on the AnalysisManager
240
241   if(!fMCEventH)return;
242   if(!fMCEventH->MCEvent())return;
243   AliStack *pStack = fMCEventH->MCEvent()->Stack();
244   if(!pStack)return;
245
246   fMCEventH->CreateLabelMap();
247
248   //
249   // Get the Event Header 
250   // 
251
252   AliHeader* header = fMCEventH->MCEvent()->Header();
253    // get the MC vertex
254   AliGenEventHeader* genHeader = 0;
255   if (header) genHeader = header->GenEventHeader();
256   if (genHeader) {
257       TArrayF vtxMC(3);
258       genHeader->PrimaryVertex(vtxMC);
259       mcHeader->SetVertex(vtxMC[0],vtxMC[1],vtxMC[2]);
260
261       // we search the MCEventHeaders first 
262       // Two cases, cocktail or not...
263       AliGenCocktailEventHeader* genCocktailHeader = dynamic_cast<AliGenCocktailEventHeader*>(genHeader);
264       if(genCocktailHeader){
265           // we have a coktail header
266           mcHeader->AddGeneratorName(genHeader->GetName());
267           // Loop from the back so that the first one sets the process type
268           TList* headerList = genCocktailHeader->GetHeaders();
269           for(int i = headerList->GetEntries()-1;i>=0;--i){
270               AliGenEventHeader *headerEntry = dynamic_cast<AliGenEventHeader*>(headerList->At(i));
271               SetMCHeaderInfo(mcHeader,headerEntry);
272           }
273       }
274       else{
275           // No Cocktail just take the first one
276           SetMCHeaderInfo(mcHeader,genHeader);
277       }
278   }
279   
280
281
282
283
284   // Store the AliAODParticlesMC
285   AliMCEvent* mcEvent = fMCEventH->MCEvent();
286   
287   Int_t np    = mcEvent->GetNumberOfTracks();
288   Int_t nprim = mcEvent->GetNumberOfPrimaries();
289
290
291   Int_t j = 0;
292   TClonesArray& l = *mcarray;
293
294   for(int i = 0; i < np; ++i){
295       if(fMCEventH->IsParticleSelected(i)){
296           Int_t flag = 0;
297           AliMCParticle* mcpart =  (AliMCParticle*) mcEvent->GetTrack(i);
298           if(i<nprim)flag |= AliAODMCParticle::kPrimary;
299           
300           if(mcEvent->IsPhysicalPrimary(i))flag |= AliAODMCParticle::kPhysicalPrim;
301           
302           if(fMCEventH->GetNewLabel(i)!=j){
303               AliError(Form("MISMATCH New label %d j: %d",fMCEventH->GetNewLabel(i),j));
304           }
305
306           AliAODMCParticle mcpart_tmp(mcpart,i,flag);
307           
308           mcpart_tmp.SetStatus(mcpart->Particle()->GetStatusCode());
309           // 
310           Int_t d0 =  mcpart_tmp.GetDaughter(0);
311           Int_t d1 =  mcpart_tmp.GetDaughter(1);
312           Int_t m =   mcpart_tmp.GetMother();
313           
314           // other than for the track labels, negative values mean
315           // no daughter/mother so preserve it
316           
317           if(d0<0 && d1<0){
318               // no first daughter -> no second daughter
319               // nothing to be done
320               // second condition not needed just for sanity check at the end
321               mcpart_tmp.SetDaughter(0,d0);
322               mcpart_tmp.SetDaughter(1,d1);
323           } else if(d1 < 0 && d0 >= 0) {
324               // Only one daughter
325               // second condition not needed just for sanity check at the end
326               if(fMCEventH->IsParticleSelected(d0)){
327                   mcpart_tmp.SetDaughter(0,fMCEventH->GetNewLabel(d0));
328               } else {
329                   mcpart_tmp.SetDaughter(0,-1);
330               }
331               mcpart_tmp.SetDaughter(1,d1);
332           }
333           else if (d0 > 0 && d1 > 0 ){
334               // we have two or more daughters loop on the stack to see if they are
335               // selected
336               Int_t d0_tmp = -1;
337               Int_t d1_tmp = -1;
338               for(int id = d0; id<=d1;++id){
339                   if(fMCEventH->IsParticleSelected(id)){
340                       if(d0_tmp==-1){
341                           // first time
342                           d0_tmp = fMCEventH->GetNewLabel(id);
343                           d1_tmp = d0_tmp; // this is to have the same schema as on the stack i.e. with one daugther d0 and d1 are the same 
344                       }
345                       else d1_tmp = fMCEventH->GetNewLabel(id);
346                   }
347               }
348               mcpart_tmp.SetDaughter(0,d0_tmp);
349               mcpart_tmp.SetDaughter(1,d1_tmp);
350           } else {
351               AliError(Form("Unxpected indices %d %d",d0,d1));
352           }
353           
354           if(m<0){
355               mcpart_tmp.SetMother(m);
356           } else {
357               if(fMCEventH->IsParticleSelected(m))mcpart_tmp.SetMother(fMCEventH->GetNewLabel(m));
358               else AliError(Form("PROBLEM Mother not selected %d \n", m));
359           }
360
361           new (l[j++]) AliAODMCParticle(mcpart_tmp);
362           
363       }
364   }
365   AliInfo(Form("AliAODHandler::StoreMCParticles: Selected %d (Primaries %d / total %d) after validation \n",
366                j,nprim,np));
367   
368   // Set the labels in the AOD output...
369   // Remapping
370
371   // AODTracks
372   TClonesArray* tracks = fAODEvent->GetTracks();
373   if(tracks){
374     for(int it = 0; it < fAODEvent->GetNTracks();++it){
375       AliAODTrack *track = fAODEvent->GetTrack(it);
376       
377       Int_t sign = 1;
378       Int_t label = track->GetLabel();
379       if(label<0){ // preserve the sign for later usage
380         label *= -1;
381         sign  = -1;
382       }
383
384       if (label >= AliMCEvent::BgLabelOffset()) label =  mcEvent->BgLabelToIndex(label);
385       if(label > np || track->GetLabel() == 0){
386         AliWarning(Form("Wrong ESD track label %5d (%5d)",track->GetLabel(), label));
387       }
388       if(fMCEventH->GetNewLabel(label) == 0){
389         AliWarning(Form("New label not found for %5d (%5d)",track->GetLabel(), label));
390       }
391       track->SetLabel(sign*fMCEventH->GetNewLabel(label));
392     }
393   }
394   
395   // AOD calo cluster
396   TClonesArray *clusters = fAODEvent->GetCaloClusters();
397   if(clusters){
398     for (Int_t iClust = 0;iClust < fAODEvent->GetNumberOfCaloClusters(); ++iClust) {
399       AliAODCaloCluster * cluster = fAODEvent->GetCaloCluster(iClust);
400       UInt_t nLabel    = cluster->GetNLabels();
401       // Ugly but do not want to fragment memory by creating 
402       // new Int_t (nLabel)
403       Int_t* labels    = const_cast<Int_t*>(cluster->GetLabels());
404       if (labels){
405         for(UInt_t i = 0;i < nLabel;++i){
406           labels[i] = fMCEventH->GetNewLabel(cluster->GetLabelAt(i));
407         }
408       }
409       //      cluster->SetLabels(labels,nLabel);
410     }// iClust
411   }// clusters
412
413   // AOD tracklets
414   AliAODTracklets *tracklets = fAODEvent->GetTracklets();
415   if(tracklets){
416     for(int it = 0;it < tracklets->GetNumberOfTracklets();++it){
417       int label0 = tracklets->GetLabel(it,0);
418       int label1 = tracklets->GetLabel(it,1);
419       if(label0>=0)label0 = fMCEventH->GetNewLabel(label0);      
420       if(label1>=0)label1 = fMCEventH->GetNewLabel(label1);
421       tracklets->SetLabel(it,0,label0);
422       tracklets->SetLabel(it,1,label1);
423     }
424   }
425
426 }
427
428 //______________________________________________________________________________
429 Bool_t AliAODHandler::FinishEvent()
430 {
431   // Fill data structures
432   
433   if(fFillAOD && fFillAODRun && fAODEvent){
434       fAODEvent->MakeEntriesReferencable();
435       fTreeA->BranchRef();
436       FillTree();
437   }
438
439   if (fFillAOD && fFillAODRun) {      
440     if (fExtensions) {
441       TIter next(fExtensions);
442       AliAODExtension *ext;
443       while ((ext=(AliAODExtension*)next())) ext->FinishEvent();
444     }
445     if (fFilters) {   
446       TIter nextf(fFilters);
447       AliAODExtension *ext;
448       while ((ext=(AliAODExtension*)nextf())) {
449               ext->FinishEvent();
450       }  
451     }       
452   }  
453   
454   if (fIsStandard) 
455   {
456     fAODEvent->ResetStd();    
457   }
458   
459   if (fAODEvent) 
460   {
461     TClonesArray *mcarray = static_cast<TClonesArray*>(fAODEvent->FindListObject(AliAODMCParticle::StdBranchName()));
462     if(mcarray) mcarray->Delete();
463     
464     AliAODMCHeader *mcHeader = static_cast<AliAODMCHeader*>(fAODEvent->FindListObject(AliAODMCHeader::StdBranchName()));
465     if(mcHeader) mcHeader->Reset();
466   }
467   
468   // Reset AOD replication flag
469   fAODIsReplicated = kFALSE;
470   return kTRUE;
471 }
472
473 //______________________________________________________________________________
474 Bool_t AliAODHandler::Terminate()
475 {
476   // Terminate 
477   AddAODtoTreeUserInfo();
478   
479   TIter nextF(fFilters);
480   AliAODExtension *ext;
481   while ((ext=static_cast<AliAODExtension*>(nextF())))
482   {
483     ext->AddAODtoTreeUserInfo();
484   }
485
486   TIter nextE(fExtensions);
487   while ((ext=static_cast<AliAODExtension*>(nextE())))
488   {
489     ext->AddAODtoTreeUserInfo();
490   }
491   
492   return kTRUE;
493 }
494
495 //______________________________________________________________________________
496 Bool_t AliAODHandler::TerminateIO()
497 {
498   // Terminate IO
499   if (fFileA) {
500     fFileA->Write();
501     fFileA->Close();
502     delete fFileA;
503     fFileA = 0;
504     // When closing the file, the tree is also deleted.
505     fTreeA = 0;
506   }
507   
508   TIter nextF(fFilters);
509   AliAODExtension *ext;
510   while ((ext=static_cast<AliAODExtension*>(nextF())))
511   {
512     ext->TerminateIO();
513   }  
514
515   TIter nextE(fExtensions);
516   while ((ext=static_cast<AliAODExtension*>(nextE())))
517   {
518     ext->TerminateIO();
519   }  
520   
521   return kTRUE;
522 }
523
524 //______________________________________________________________________________
525 void AliAODHandler::CreateTree(Int_t flag)
526 {
527     // Creates the AOD Tree
528     fTreeA = new TTree("aodTree", "AliAOD tree");
529     fTreeA->Branch(fAODEvent->GetList());
530     if (flag == 0) fTreeA->SetDirectory(0);
531 }
532
533 //______________________________________________________________________________
534 void AliAODHandler::FillTree()
535 {
536  
537     // Fill the AOD Tree
538     fTreeA->Fill();
539 }
540
541 //______________________________________________________________________________
542 void AliAODHandler::AddAODtoTreeUserInfo()
543 {
544   // Add aod event to tree user info
545   if (fTreeA) fTreeA->GetUserInfo()->Add(fAODEvent);
546   // Now the tree owns our fAODEvent...
547   fAODEvent = 0;
548 }
549
550 //______________________________________________________________________________
551 void AliAODHandler::AddBranch(const char* cname, void* addobj, const char* filename)
552 {
553   // Add a new branch to the aod. Added optional filename parameter if the
554   // branch should be written to a separate file.
555   
556   if (strlen(filename)) 
557   {
558     AliAODExtension *ext = AddExtension(filename);
559     ext->AddBranch(cname, addobj);
560     return;
561   } 
562   
563   // Add branch to all filters
564   // Add branch to all filters
565   if (fFilters) {
566     TIter next(fFilters);
567     AliAODExtension *ext;
568     while ((ext=(AliAODExtension*)next())) ext->AddBranch(cname, addobj);
569   }
570   
571   TDirectory *owd = gDirectory;
572   if (fFileA) 
573   {
574     fFileA->cd();
575   }
576
577   char** apointer = (char**) addobj;
578   TObject* obj = (TObject*) *apointer;
579   
580   fAODEvent->AddObject(obj);
581   
582   const Int_t kSplitlevel = 99; // default value in TTree::Branch()
583   const Int_t kBufsize = 32000; // default value in TTree::Branch()
584   
585   if (!fTreeA->FindBranch(obj->GetName())) 
586   {
587     // Do the same as if we book via 
588     // TTree::Branch(TCollection*)
589     
590     fTreeA->Bronch(obj->GetName(), cname, fAODEvent->GetList()->GetObjectRef(obj),
591                    kBufsize, kSplitlevel - 1);
592   }
593   owd->cd();
594 }
595
596 //______________________________________________________________________________
597 AliAODExtension *AliAODHandler::AddExtension(const char *filename, const char *title)
598 {
599   // Add an AOD extension with some branches in a different file.
600   
601   TString fname(filename);
602   if (!fname.EndsWith(".root")) fname += ".root";
603   if (!fExtensions) {
604     fExtensions = new TObjArray();
605     fExtensions->SetOwner();
606   }   
607   AliAODExtension *ext = (AliAODExtension*)fExtensions->FindObject(fname);
608   if (!ext) {
609     ext = new AliAODExtension(fname, title);
610     fExtensions->Add(ext);
611   }   
612   return ext;
613 }
614
615 //______________________________________________________________________________
616 AliAODExtension *AliAODHandler::GetExtension(const char *filename) const
617 {
618   // Getter for AOD extensions via file name.
619   if (!fExtensions) return NULL;
620   return (AliAODExtension*)fExtensions->FindObject(filename);
621 }   
622
623 //______________________________________________________________________________
624 AliAODExtension *AliAODHandler::AddFilteredAOD(const char *filename, const char *filtername)
625 {
626   // Add an AOD extension that can write only AOD events that pass a user filter.
627   if (!fFilters) {
628     fFilters = new TObjArray();
629     fFilters->SetOwner();
630   } 
631   AliAODExtension *filter = (AliAODExtension*)fFilters->FindObject(filename);
632   if (!filter) {
633     filter = new AliAODExtension(filename, filtername, kTRUE);
634     fFilters->Add(filter);
635   }
636   return filter;
637 }      
638
639 //______________________________________________________________________________
640 AliAODExtension *AliAODHandler::GetFilteredAOD(const char *filename) const
641 {
642   // Getter for AOD filters via file name.
643   if (!fFilters) return NULL;
644   return (AliAODExtension*)fFilters->FindObject(filename);
645 }   
646
647 //______________________________________________________________________________
648 void AliAODHandler::SetOutputFileName(const char* fname)
649 {
650 // Set file name.
651    fFileName = fname;
652 }
653
654 //______________________________________________________________________________
655 const char *AliAODHandler::GetOutputFileName()
656 {
657 // Get file name.
658    return fFileName.Data();
659 }
660
661 //______________________________________________________________________________
662 const char *AliAODHandler::GetExtraOutputs() const
663 {
664   // Get extra outputs as a string separated by commas.
665   static TString eoutputs;
666   eoutputs = "";
667   TObject *obj;
668   if (fExtensions) {
669     TIter next1(fExtensions);
670     while ((obj=next1())) {
671       if (!eoutputs.IsNull()) eoutputs += ",";
672       eoutputs += obj->GetName();
673     }
674   }
675   if (fFilters) {
676     TIter next2(fFilters);
677     while ((obj=next2())) {
678       if (!eoutputs.IsNull()) eoutputs += ",";
679       eoutputs += obj->GetName();
680     }
681   }
682   return eoutputs.Data();
683 }
684
685 //______________________________________________________________________________
686 Bool_t AliAODHandler::HasExtensions() const
687 {
688   // Whether or not we manage extensions
689   
690   if ( fExtensions && fExtensions->GetEntries()>0 ) return kTRUE;
691   
692   return kFALSE;
693 }
694
695 //______________________________________________________________________________
696 void  AliAODHandler::SetMCHeaderInfo(AliAODMCHeader *mcHeader,AliGenEventHeader *genHeader){
697
698
699   // Utility function to cover different cases for the AliGenEventHeader
700   // Needed since different ProcessType and ImpactParamter are not 
701   // in the base class...
702   // We don't encode process types for event cocktails yet
703   // could be done e.g. by adding offsets depnding on the generator
704
705   mcHeader->AddGeneratorName(genHeader->GetName());
706   if(!genHeader)return;
707   AliGenPythiaEventHeader *pythiaGenHeader = dynamic_cast<AliGenPythiaEventHeader*>(genHeader);
708     if (pythiaGenHeader) {
709       mcHeader->SetEventType(pythiaGenHeader->ProcessType());
710       mcHeader->SetPtHard(pythiaGenHeader->GetPtHard());
711       return;
712     }
713     
714     AliGenDPMjetEventHeader* dpmJetGenHeader = dynamic_cast<AliGenDPMjetEventHeader*>(genHeader);
715
716   if (dpmJetGenHeader){
717     mcHeader->SetEventType(dpmJetGenHeader->ProcessType());
718     return;
719   } 
720
721   AliGenHijingEventHeader* hijingGenHeader = dynamic_cast<AliGenHijingEventHeader*>(genHeader);
722   if(hijingGenHeader){
723     mcHeader->SetImpactParameter(hijingGenHeader->ImpactParameter());
724     return;
725   }
726   
727   AliWarning(Form("MC Eventheader not known: %s",genHeader->GetName()));
728
729 }
730