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