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