]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliAODHandler.cxx
Fixed initialization of AOD file in proof case. Handling delta AOD's supported.
[u/mrichter/AliRoot.git] / STEER / AliAODHandler.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2007, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 //-------------------------------------------------------------------------
19 //     Implementation of the Virtual Event Handler Interface for AOD
20 //     Author: Andreas Morsch, CERN
21 //-------------------------------------------------------------------------
22
23
24 #include <TTree.h>
25 #include <TFile.h>
26 #include <TRef.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     fAODIsReplicated(kFALSE),
65     fAODEvent(NULL),
66     fMCEventH(NULL),
67     fTreeA(NULL),
68     fFileA(NULL),
69     fFileName(""),
70     fExtensions(NULL)
71 {
72   // default constructor
73 }
74
75 //______________________________________________________________________________
76 AliAODHandler::AliAODHandler(const char* name, const char* title):
77     AliVEventHandler(name, title),
78     fIsStandard(kTRUE),
79     fFillAOD(kTRUE),
80     fNeedsHeaderReplication(kFALSE),
81     fNeedsTracksBranchReplication(kFALSE),
82     fNeedsVerticesBranchReplication(kFALSE),
83     fNeedsV0sBranchReplication(kFALSE),
84     fNeedsTrackletsBranchReplication(kFALSE),
85     fNeedsPMDClustersBranchReplication(kFALSE),
86     fNeedsJetsBranchReplication(kFALSE),
87     fNeedsFMDClustersBranchReplication(kFALSE),
88     fNeedsCaloClustersBranchReplication(kFALSE),
89     fAODIsReplicated(kFALSE),
90     fAODEvent(NULL),
91     fMCEventH(NULL),
92     fTreeA(NULL),
93     fFileA(NULL),
94     fFileName(""),
95     fExtensions(NULL)
96 {
97 }
98
99 //______________________________________________________________________________
100 AliAODHandler::~AliAODHandler() 
101 {
102  // Destructor.
103   delete fAODEvent;
104   if(fFileA){
105     // is already handled in TerminateIO
106     fFileA->Close();
107     delete fFileA;
108   }
109   delete fTreeA;
110   if (fExtensions) delete fExtensions;
111 }
112
113 //______________________________________________________________________________
114 Bool_t AliAODHandler::Init(Option_t* opt)
115 {
116   // Initialize IO
117   //
118   // Create the AODevent object
119   if(!fAODEvent){
120     fAODEvent = new AliAODEvent();
121     if (fIsStandard) fAODEvent->CreateStdContent();
122   }
123   //
124   // File opening according to execution mode
125   TString option(opt);
126   option.ToLower();
127   TDirectory *owd = gDirectory;
128   if (option.Contains("proof")) {
129     // proof
130     // Merging via files. Need to access analysis manager via interpreter.
131     gROOT->ProcessLine(Form("AliAnalysisManager::GetAnalysisManager()->OpenProofFile(%s, \"RECREATE\");", fFileName.Data()));
132     gROOT->ProcessLine(Form("AliAnalysisManager::GetAnalysisManager()->GetCommonOutputContainer()->SetFile((TFile*)0x%lx);", gFile));
133     fFileA = gFile;
134   } else {
135     // local and grid
136     fFileA = new TFile(fFileName.Data(), "RECREATE");
137   }
138   CreateTree(1);
139   owd->cd();
140   if (fExtensions) {
141      TIter next(fExtensions);
142      AliAODExtension *ext;
143      while ((ext=(AliAODExtension*)next())) ext->Init(option);
144   }   
145   return kTRUE;
146 }
147
148 //______________________________________________________________________________
149 void AliAODHandler::StoreMCParticles(){
150
151   // 
152   // Remap the labels from ESD stack and store
153   // the AODMCParticles, makes only sense if we have
154   // the mcparticles branch
155   // has to be done here since we cannot know in advance 
156   // which particles are needed (e.g. by the tracks etc.)
157   //
158   // Particles have been selected by AliMCEventhanlder->SelectParticle()
159   // To use the MCEventhandler here we need to set it from the outside
160   // can vanish when Handler go to the ANALYSISalice library
161   //
162   // The Branch booking for mcParticles and mcHeader has to happen 
163   // in an external task for now since the AODHandler does not have access
164   // the AnalysisManager. For the same reason the pointer t o the MCEventH
165   // has to passed to the AOD Handler by this task 
166   // (doing this in the steering macro would not work on PROOF)
167
168   TClonesArray *mcarray = (TClonesArray*)fAODEvent->FindListObject(AliAODMCParticle::StdBranchName()); 
169   if(!mcarray)return;
170   mcarray->Delete();
171
172   AliAODMCHeader *mcHeader = (AliAODMCHeader*)fAODEvent->FindListObject(AliAODMCHeader::StdBranchName()); 
173   if(!mcHeader)return;
174
175   mcHeader->Reset();
176
177   // Get the MC Infos.. Handler needs to be set before 
178   // while adding the branch
179   // This needs to be done, not to depend on the AnalysisManager
180
181   if(!fMCEventH)return;
182   if(!fMCEventH->MCEvent())return;
183   AliStack *pStack = fMCEventH->MCEvent()->Stack();
184   if(!pStack)return;
185
186   fMCEventH->CreateLabelMap();
187
188   //
189   // Get the Event Header 
190   // 
191
192   AliHeader* header = fMCEventH->MCEvent()->Header();
193   if (!header)return;
194
195   // get the MC vertex
196   AliGenEventHeader* genHeader = header->GenEventHeader();
197   TArrayF vtxMC(3);
198   genHeader->PrimaryVertex(vtxMC);
199   mcHeader->SetVertex(vtxMC[0],vtxMC[1],vtxMC[2]);
200
201   // we search the MCEventHeaders first 
202   // Two cases, cocktail or not...
203   AliGenCocktailEventHeader* genCocktailHeader = dynamic_cast<AliGenCocktailEventHeader*>(genHeader);
204   if(genCocktailHeader){
205     // we have a coktail header
206     mcHeader->AddGeneratorName(genHeader->GetName());
207     // Loop from the back so that the first one sets the process type
208     TList* headerList = genCocktailHeader->GetHeaders();
209     for(int i = headerList->GetEntries()-1;i>=0;--i){
210       AliGenEventHeader *headerEntry = dynamic_cast<AliGenEventHeader*>(headerList->At(i));
211       SetMCHeaderInfo(mcHeader,headerEntry);
212     }
213   }
214   else{
215     // No Cocktail just take the first one
216     SetMCHeaderInfo(mcHeader,genHeader);
217   }
218
219
220
221
222   // Store the AliAODParticlesMC
223
224   Int_t np    = pStack->GetNtrack();
225   Int_t nprim = pStack->GetNprimary();
226
227
228   Int_t j = 0;
229   TClonesArray& l = *mcarray;
230
231   for(int i = 0;i < np;++i){
232     if(fMCEventH->IsParticleSelected(i)){
233
234       Int_t flag = 0;
235       TParticle *part = pStack->Particle(i);
236       if(i<nprim)flag |= AliAODMCParticle::kPrimary;
237       if(pStack->IsPhysicalPrimary(i))flag |= AliAODMCParticle::kPhysicalPrim;
238
239       if(fMCEventH->GetNewLabel(i)!=j){
240         AliError(Form("MISMATCH New label %d j: %d",fMCEventH->GetNewLabel(i),j));
241       }
242       AliAODMCParticle mcpart_tmp(part,i,flag);
243
244       // 
245       Int_t d0 =  mcpart_tmp.GetDaughter(0);
246       Int_t d1 =  mcpart_tmp.GetDaughter(1);
247       Int_t m =  mcpart_tmp.GetMother();
248
249       // other than for the track labels, negative values mean
250       // no daughter/mother so preserve it
251  
252       if(d0<0 && d1<0){
253         // no first daughter -> no second daughter
254         // nothing to be done
255         // second condition not needed just for sanity check at the end
256         mcpart_tmp.SetDaughter(0,d0);
257         mcpart_tmp.SetDaughter(1,d1);
258       }
259       else if(d1 < 0 && d0 >= 0){
260         // Only one daughter
261         // second condition not needed just for sanity check at the end
262         if(fMCEventH->IsParticleSelected(d0)){
263           mcpart_tmp.SetDaughter(0,fMCEventH->GetNewLabel(d0));
264         }
265         else{
266           mcpart_tmp.SetDaughter(0,-1);
267         }
268         mcpart_tmp.SetDaughter(1,d1);
269       }
270       else if (d0 > 0 && d1 > 0 ){
271         // we have two or more daughters loop on the stack to see if they are
272         // selected
273         Int_t d0_tmp = -1;
274         Int_t d1_tmp = -1;
275         for(int id = d0; id<=d1;++id){
276           if(fMCEventH->IsParticleSelected(id)){
277             if(d0_tmp==-1){
278               // first time
279               d0_tmp = fMCEventH->GetNewLabel(id);
280               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 
281             }
282             else d1_tmp = fMCEventH->GetNewLabel(id);
283           }
284         }
285         mcpart_tmp.SetDaughter(0,d0_tmp);
286         mcpart_tmp.SetDaughter(1,d1_tmp);
287       }
288       else{
289         AliError(Form("Unxpected indices %d %d",d0,d1));
290       }
291
292       if(m<0){
293         mcpart_tmp.SetMother(m);
294       }
295       else{
296         if(fMCEventH->IsParticleSelected(m))mcpart_tmp.SetMother(fMCEventH->GetNewLabel(m));
297         else AliError("PROBLEM Mother not selected");
298       }
299
300       new (l[j++]) AliAODMCParticle(mcpart_tmp);
301       
302     }
303   }
304   AliInfo(Form("AliAODHandler::StoreMCParticles: Selected %d (Primaries %d / total %d) after validation",
305                j,nprim,np));
306   
307   // Set the labels in the AOD output...
308   // Remapping
309
310   // AODTracks
311   TClonesArray* tracks = fAODEvent->GetTracks();
312   if(tracks){
313     for(int it = 0; it < fAODEvent->GetNTracks();++it){
314       AliAODTrack *track = fAODEvent->GetTrack(it);
315       
316       if(TMath::Abs(track->GetLabel())>np||track->GetLabel()==0){
317         AliWarning(Form("Wrong ESD track label %d",track->GetLabel()));
318       }
319       if(fMCEventH->GetNewLabel(track->GetLabel())==0){
320         AliWarning(Form("New label not found for %d",track->GetLabel()));
321       }
322       track->SetLabel(fMCEventH->GetNewLabel(track->GetLabel()));
323     }
324   }
325   
326   // AOD calo cluster
327   TClonesArray *clusters = fAODEvent->GetCaloClusters();
328   if(clusters){
329     for (Int_t iClust = 0;iClust < fAODEvent->GetNCaloClusters(); ++iClust) {
330       AliAODCaloCluster * cluster = fAODEvent->GetCaloCluster(iClust);
331       UInt_t nLabel    = cluster->GetNLabel();
332       // Ugly but do not want to fragment memory by creating 
333       // new Int_t (nLabel)
334       Int_t* labels    = const_cast<Int_t*>(cluster->GetLabels());
335       if (labels){
336         for(UInt_t i = 0;i < nLabel;++i){
337           labels[i] = fMCEventH->GetNewLabel(cluster->GetLabel(i));
338         }
339       }
340       //      cluster->SetLabels(labels,nLabel);
341     }// iClust
342   }// clusters
343
344   // AOD tracklets
345   AliAODTracklets *tracklets = fAODEvent->GetTracklets();
346   if(tracklets){
347     for(int it = 0;it < tracklets->GetNumberOfTracklets();++it){
348       int label0 = tracklets->GetLabel(it,0);
349       int label1 = tracklets->GetLabel(it,1);
350       if(label0>=0)label0 = fMCEventH->GetNewLabel(label0);      
351       if(label1>=0)label1 = fMCEventH->GetNewLabel(label1);
352       tracklets->SetLabel(it,0,label0);
353       tracklets->SetLabel(it,1,label1);
354     }
355   }
356
357 }
358
359 //______________________________________________________________________________
360 Bool_t AliAODHandler::FinishEvent()
361 {
362   // Fill data structures
363   if(fFillAOD){
364     fAODEvent->MakeEntriesReferencable();
365     StoreMCParticles();
366     FillTree();
367     if (fExtensions) {
368       TIter next(fExtensions);
369       AliAODExtension *ext;
370       while ((ext=(AliAODExtension*)next())) {
371         ext->GetAOD()->MakeEntriesReferencable();
372         ext->GetTree()->Fill();
373       }  
374     }       
375   }
376
377   if (fIsStandard) fAODEvent->ResetStd();
378   // Reset AOD replication flag
379   fAODIsReplicated = kFALSE;
380   return kTRUE;
381 }
382
383 //______________________________________________________________________________
384 Bool_t AliAODHandler::Terminate()
385 {
386   // Terminate 
387   AddAODtoTreeUserInfo();
388   if (fExtensions) {
389     TIter next(fExtensions);
390     AliAODExtension *ext;
391     while ((ext=(AliAODExtension*)next())) ext->GetTree()->GetUserInfo()->Add(ext->GetAOD());
392   }  
393   return kTRUE;
394 }
395
396 //______________________________________________________________________________
397 Bool_t AliAODHandler::TerminateIO()
398 {
399   // Terminate IO
400   if (fFileA) {
401     fFileA->Close();
402     delete fFileA;
403     fFileA = 0;
404   }
405   if (fExtensions) {
406     TIter next(fExtensions);
407     AliAODExtension *ext;
408     while ((ext=(AliAODExtension*)next())) ext->TerminateIO();
409   }  
410   return kTRUE;
411 }
412
413 //______________________________________________________________________________
414 void AliAODHandler::CreateTree(Int_t flag)
415 {
416     // Creates the AOD Tree
417     fTreeA = new TTree("aodTree", "AliAOD tree");
418     fTreeA->Branch(fAODEvent->GetList());
419     TRef junk = (TObject*)fTreeA->BranchRef();
420     if (flag == 0) fTreeA->SetDirectory(0);
421 }
422
423 //______________________________________________________________________________
424 void AliAODHandler::FillTree()
425 {
426     // Fill the AOD Tree
427   fTreeA->Fill();
428 }
429
430 //______________________________________________________________________________
431 void AliAODHandler::AddAODtoTreeUserInfo()
432 {
433   // Add aod event to tree user info
434   fTreeA->GetUserInfo()->Add(fAODEvent);
435 }
436
437 //______________________________________________________________________________
438 void AliAODHandler::AddBranch(const char* cname, void* addobj, const char* filename)
439 {
440     // Add a new branch to the aod. Added optional filename parameter if the
441     // branch should be written to a separate file.
442     if (strlen(filename)) {
443        AliAODExtension *ext = AddExtension(filename);
444        ext->AddBranch(cname, addobj);
445        return;
446     }
447     TDirectory *owd = gDirectory;
448     if (fFileA) {
449       fFileA->cd();
450     }
451     char** apointer = (char**) addobj;
452     TObject* obj = (TObject*) *apointer;
453
454     fAODEvent->AddObject(obj);
455  
456     const Int_t kSplitlevel = 99; // default value in TTree::Branch()
457     const Int_t kBufsize = 32000; // default value in TTree::Branch()
458
459     if (!fTreeA->FindBranch(obj->GetName())) {
460       // Do the same as if we book via 
461       // TTree::Branch(TCollection*)
462       
463       fTreeA->Bronch(obj->GetName(), cname, fAODEvent->GetList()->GetObjectRef(obj),
464                      kBufsize, kSplitlevel - 1);
465       //    fTreeA->Branch(obj->GetName(), cname, addobj);
466     }
467     owd->cd();
468 }
469
470 //______________________________________________________________________________
471 AliAODExtension *AliAODHandler::AddExtension(const char *filename, const char *title)
472 {
473 // Add an AOD extension with some branches in a different file.
474    TString fname(filename);
475    if (!fname.EndsWith(".root")) fname += ".root";
476    if (!fExtensions) {
477       fExtensions = new TObjArray();
478       fExtensions->SetOwner();
479    }   
480    AliAODExtension *ext = (AliAODExtension*)fExtensions->FindObject(fname);
481    if (!ext) {
482       ext = new AliAODExtension(fname, title);
483       fExtensions->Add(ext);
484    }   
485    return ext;
486 }
487    
488 //______________________________________________________________________________
489 void AliAODHandler::SetOutputFileName(const char* fname)
490 {
491 // Set file name.
492    fFileName = fname;
493 }
494
495 //______________________________________________________________________________
496 const char *AliAODHandler::GetOutputFileName()
497 {
498 // Get file name.
499    return fFileName.Data();
500 }
501
502 //______________________________________________________________________________
503 void  AliAODHandler::SetMCHeaderInfo(AliAODMCHeader *mcHeader,AliGenEventHeader *genHeader){
504
505
506   // Utility function to cover different cases for the AliGenEventHeader
507   // Needed since different ProcessType and ImpactParamter are not 
508   // in the base class...
509   // We don't encode process types for event cocktails yet
510   // could be done e.g. by adding offsets depnding on the generator
511
512   mcHeader->AddGeneratorName(genHeader->GetName());
513   if(!genHeader)return;
514   AliGenPythiaEventHeader *pythiaGenHeader = dynamic_cast<AliGenPythiaEventHeader*>(genHeader);
515     if (pythiaGenHeader) {
516       mcHeader->SetEventType(pythiaGenHeader->ProcessType());
517       return;
518     }
519     
520     AliGenDPMjetEventHeader* dpmJetGenHeader = dynamic_cast<AliGenDPMjetEventHeader*>(genHeader);
521
522   if (dpmJetGenHeader){
523     mcHeader->SetEventType(dpmJetGenHeader->ProcessType());
524     return;
525   } 
526
527   AliGenHijingEventHeader* hijingGenHeader = dynamic_cast<AliGenHijingEventHeader*>(genHeader);
528   if(hijingGenHeader){
529     mcHeader->SetImpactParameter(hijingGenHeader->ImpactParameter());
530     return;
531   }
532   
533   AliWarning(Form("MC Eventheader not known: %s",genHeader->GetName()));
534
535 }
536
537 ClassImp(AliAODExtension)
538
539 //-------------------------------------------------------------------------
540 //     Support class for AOD extensions. This is created by the user analysis
541 //     that requires a separate file for some AOD branches. The name of the 
542 //     AliAODExtension object is the file name where the AOD branches will be
543 //     stored.
544 //-------------------------------------------------------------------------
545
546 //______________________________________________________________________________
547 AliAODExtension::~AliAODExtension()
548 {
549 // Destructor.
550   delete fAODEvent;
551   if(fFileE){
552     // is already handled in TerminateIO
553     fFileE->Close();
554     delete fFileE;
555   }
556   delete fTreeE;
557 }
558
559 //______________________________________________________________________________
560 void AliAODExtension::AddBranch(const char* cname, void* addobj)
561 {
562     // Add a new branch to the aod 
563     if (!fAODEvent) Init("");
564     TDirectory *owd = gDirectory;
565     if (fFileE) {
566       fFileE->cd();
567     }
568     char** apointer = (char**) addobj;
569     TObject* obj = (TObject*) *apointer;
570
571     fAODEvent->AddObject(obj);
572  
573     const Int_t kSplitlevel = 99; // default value in TTree::Branch()
574     const Int_t kBufsize = 32000; // default value in TTree::Branch()
575
576     if (!fTreeE->FindBranch(obj->GetName())) {
577       // Do the same as if we book via 
578       // TTree::Branch(TCollection*)
579       
580       fTreeE->Bronch(obj->GetName(), cname, fAODEvent->GetList()->GetObjectRef(obj),
581                      kBufsize, kSplitlevel - 1);
582       //    fTreeA->Branch(obj->GetName(), cname, addobj);
583     }
584     owd->cd();
585 }
586
587 //______________________________________________________________________________
588 Bool_t AliAODExtension::Init(Option_t *option)
589 {
590 // Initialize IO.
591   if(!fAODEvent) fAODEvent = new AliAODEvent();
592   TDirectory *owd = gDirectory;
593   TString opt(option);
594   opt.ToLower();
595   if (opt.Contains("proof")) {
596     // proof
597     // Merging via files. Need to access analysis manager via interpreter.
598     gROOT->ProcessLine(Form("AliAnalysisManager::GetAnalysisManager()->OpenProofFile(%s, \"RECREATE\");", fName.Data()));
599     fFileE = gFile;
600   } else {
601     fFileE = new TFile(GetName(), "RECREATE");
602   }  
603   fTreeE = new TTree("aodTree", "AliAOD tree");
604   fTreeE->Branch(fAODEvent->GetList());
605   TRef junk = (TObject*)fTreeE->BranchRef();
606   owd->cd();
607   return kTRUE;
608 }
609
610 //______________________________________________________________________________
611 Bool_t AliAODExtension::TerminateIO()
612 {
613   // Terminate IO
614   if (fFileE) {
615     fFileE->Write();
616     fFileE->Close();
617     delete fFileE;
618     fFileE = 0;
619   }
620   return kTRUE;
621 }