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