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