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