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