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