2 /**************************************************************************
3 * Copyright(c) 1998-2007, ALICE Experiment at CERN, All rights reserved. *
5 * Author: The ALICE Off-line Project. *
6 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
19 //-------------------------------------------------------------------------
20 // Implementation of the Virtual Event Handler Interface for AOD
21 // Author: Andreas Morsch, CERN
22 //-------------------------------------------------------------------------
32 #include "AliAODHandler.h"
33 #include "AliAODEvent.h"
34 #include "AliAODTracklets.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"
48 ClassImp(AliAODHandler)
50 //______________________________________________________________________________
51 AliAODHandler::AliAODHandler() :
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),
72 // default constructor
75 //______________________________________________________________________________
76 AliAODHandler::AliAODHandler(const char* name, const char* title):
77 AliVEventHandler(name, title),
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),
99 //______________________________________________________________________________
100 AliAODHandler::~AliAODHandler()
105 // is already handled in TerminateIO
110 if (fExtensions) delete fExtensions;
113 //______________________________________________________________________________
114 Bool_t AliAODHandler::Init(Option_t* opt)
118 // Create the AODevent object
120 fAODEvent = new AliAODEvent();
121 if (fIsStandard) fAODEvent->CreateStdContent();
124 // File opening according to execution mode
127 TDirectory *owd = gDirectory;
128 if (option.Contains("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));
136 fFileA = new TFile(fFileName.Data(), "RECREATE");
141 TIter next(fExtensions);
142 AliAODExtension *ext;
143 while ((ext=(AliAODExtension*)next())) ext->Init(option);
148 //______________________________________________________________________________
149 void AliAODHandler::StoreMCParticles(){
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.)
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
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)
168 TClonesArray *mcarray = (TClonesArray*)fAODEvent->FindListObject(AliAODMCParticle::StdBranchName());
172 AliAODMCHeader *mcHeader = (AliAODMCHeader*)fAODEvent->FindListObject(AliAODMCHeader::StdBranchName());
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
181 if(!fMCEventH)return;
182 if(!fMCEventH->MCEvent())return;
183 AliStack *pStack = fMCEventH->MCEvent()->Stack();
186 fMCEventH->CreateLabelMap();
189 // Get the Event Header
192 AliHeader* header = fMCEventH->MCEvent()->Header();
196 AliGenEventHeader* genHeader = header->GenEventHeader();
198 genHeader->PrimaryVertex(vtxMC);
199 mcHeader->SetVertex(vtxMC[0],vtxMC[1],vtxMC[2]);
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);
215 // No Cocktail just take the first one
216 SetMCHeaderInfo(mcHeader,genHeader);
222 // Store the AliAODParticlesMC
223 AliMCEvent* mcEvent = fMCEventH->MCEvent();
225 Int_t np = mcEvent->GetNumberOfTracks();
226 Int_t nprim = mcEvent->GetNumberOfPrimaries();
230 TClonesArray& l = *mcarray;
232 for(int i = 0; i < np; ++i){
233 if(fMCEventH->IsParticleSelected(i)){
235 AliMCParticle* mcpart = mcEvent->GetTrack(i);
236 TParticle *part = mcpart->Particle();
237 if(i<nprim)flag |= AliAODMCParticle::kPrimary;
239 if(mcEvent->IsPhysicalPrimary(i))flag |= AliAODMCParticle::kPhysicalPrim;
241 if(fMCEventH->GetNewLabel(i)!=j){
242 AliError(Form("MISMATCH New label %d j: %d",fMCEventH->GetNewLabel(i),j));
245 AliAODMCParticle mcpart_tmp(mcpart,i,flag);
248 Int_t d0 = mcpart_tmp.GetDaughter(0);
249 Int_t d1 = mcpart_tmp.GetDaughter(1);
250 Int_t m = mcpart_tmp.GetMother();
252 // other than for the track labels, negative values mean
253 // no daughter/mother so preserve it
256 // no first daughter -> no second daughter
257 // nothing to be done
258 // second condition not needed just for sanity check at the end
259 mcpart_tmp.SetDaughter(0,d0);
260 mcpart_tmp.SetDaughter(1,d1);
261 } else if(d1 < 0 && d0 >= 0) {
263 // second condition not needed just for sanity check at the end
264 if(fMCEventH->IsParticleSelected(d0)){
265 mcpart_tmp.SetDaughter(0,fMCEventH->GetNewLabel(d0));
267 mcpart_tmp.SetDaughter(0,-1);
269 mcpart_tmp.SetDaughter(1,d1);
271 else if (d0 > 0 && d1 > 0 ){
272 // we have two or more daughters loop on the stack to see if they are
276 for(int id = d0; id<=d1;++id){
277 if(fMCEventH->IsParticleSelected(id)){
280 d0_tmp = fMCEventH->GetNewLabel(id);
281 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
283 else d1_tmp = fMCEventH->GetNewLabel(id);
286 mcpart_tmp.SetDaughter(0,d0_tmp);
287 mcpart_tmp.SetDaughter(1,d1_tmp);
289 AliError(Form("Unxpected indices %d %d",d0,d1));
293 mcpart_tmp.SetMother(m);
295 if(fMCEventH->IsParticleSelected(m))mcpart_tmp.SetMother(fMCEventH->GetNewLabel(m));
296 else AliError(Form("PROBLEM Mother not selected %d \n", m));
299 new (l[j++]) AliAODMCParticle(mcpart_tmp);
303 AliInfo(Form("AliAODHandler::StoreMCParticles: Selected %d (Primaries %d / total %d) after validation \n",
306 // Set the labels in the AOD output...
310 TClonesArray* tracks = fAODEvent->GetTracks();
312 for(int it = 0; it < fAODEvent->GetNTracks();++it){
313 AliAODTrack *track = fAODEvent->GetTrack(it);
315 Int_t label = TMath::Abs(track->GetLabel());
316 if (label >= AliMCEvent::BgLabelOffset()) label = mcEvent->BgLabelToIndex(label);
319 if(label > np || track->GetLabel() == 0){
320 AliWarning(Form("Wrong ESD track label %5d (%5d)",track->GetLabel(), label));
322 if(fMCEventH->GetNewLabel(label) == 0){
323 AliWarning(Form("New label not found for %5d (%5d)",track->GetLabel(), label));
325 track->SetLabel(fMCEventH->GetNewLabel(label));
330 TClonesArray *clusters = fAODEvent->GetCaloClusters();
332 for (Int_t iClust = 0;iClust < fAODEvent->GetNCaloClusters(); ++iClust) {
333 AliAODCaloCluster * cluster = fAODEvent->GetCaloCluster(iClust);
334 UInt_t nLabel = cluster->GetNLabel();
335 // Ugly but do not want to fragment memory by creating
336 // new Int_t (nLabel)
337 Int_t* labels = const_cast<Int_t*>(cluster->GetLabels());
339 for(UInt_t i = 0;i < nLabel;++i){
340 labels[i] = fMCEventH->GetNewLabel(cluster->GetLabel(i));
343 // cluster->SetLabels(labels,nLabel);
348 AliAODTracklets *tracklets = fAODEvent->GetTracklets();
350 for(int it = 0;it < tracklets->GetNumberOfTracklets();++it){
351 int label0 = tracklets->GetLabel(it,0);
352 int label1 = tracklets->GetLabel(it,1);
353 if(label0>=0)label0 = fMCEventH->GetNewLabel(label0);
354 if(label1>=0)label1 = fMCEventH->GetNewLabel(label1);
355 tracklets->SetLabel(it,0,label0);
356 tracklets->SetLabel(it,1,label1);
362 //______________________________________________________________________________
363 Bool_t AliAODHandler::FinishEvent()
365 // Fill data structures
367 fAODEvent->MakeEntriesReferencable();
371 TIter next(fExtensions);
372 AliAODExtension *ext;
373 while ((ext=(AliAODExtension*)next())) {
374 ext->GetAOD()->MakeEntriesReferencable();
375 ext->GetTree()->Fill();
380 if (fIsStandard) fAODEvent->ResetStd();
381 // Reset AOD replication flag
382 fAODIsReplicated = kFALSE;
386 //______________________________________________________________________________
387 Bool_t AliAODHandler::Terminate()
390 AddAODtoTreeUserInfo();
392 TIter next(fExtensions);
393 AliAODExtension *ext;
394 while ((ext=(AliAODExtension*)next())) ext->GetTree()->GetUserInfo()->Add(ext->GetAOD());
399 //______________________________________________________________________________
400 Bool_t AliAODHandler::TerminateIO()
409 TIter next(fExtensions);
410 AliAODExtension *ext;
411 while ((ext=(AliAODExtension*)next())) ext->TerminateIO();
416 //______________________________________________________________________________
417 void AliAODHandler::CreateTree(Int_t flag)
419 // Creates the AOD Tree
420 fTreeA = new TTree("aodTree", "AliAOD tree");
421 fTreeA->Branch(fAODEvent->GetList());
422 if (flag == 0) fTreeA->SetDirectory(0);
425 //______________________________________________________________________________
426 void AliAODHandler::FillTree()
432 //______________________________________________________________________________
433 void AliAODHandler::AddAODtoTreeUserInfo()
435 // Add aod event to tree user info
436 fTreeA->GetUserInfo()->Add(fAODEvent);
439 //______________________________________________________________________________
440 void AliAODHandler::AddBranch(const char* cname, void* addobj, const char* filename)
442 // Add a new branch to the aod. Added optional filename parameter if the
443 // branch should be written to a separate file.
444 if (strlen(filename)) {
445 AliAODExtension *ext = AddExtension(filename);
446 ext->AddBranch(cname, addobj);
449 TDirectory *owd = gDirectory;
453 char** apointer = (char**) addobj;
454 TObject* obj = (TObject*) *apointer;
456 fAODEvent->AddObject(obj);
458 const Int_t kSplitlevel = 99; // default value in TTree::Branch()
459 const Int_t kBufsize = 32000; // default value in TTree::Branch()
461 if (!fTreeA->FindBranch(obj->GetName())) {
462 // Do the same as if we book via
463 // TTree::Branch(TCollection*)
465 fTreeA->Bronch(obj->GetName(), cname, fAODEvent->GetList()->GetObjectRef(obj),
466 kBufsize, kSplitlevel - 1);
467 // fTreeA->Branch(obj->GetName(), cname, addobj);
472 //______________________________________________________________________________
473 AliAODExtension *AliAODHandler::AddExtension(const char *filename, const char *title)
475 // Add an AOD extension with some branches in a different file.
476 TString fname(filename);
477 if (!fname.EndsWith(".root")) fname += ".root";
479 fExtensions = new TObjArray();
480 fExtensions->SetOwner();
482 AliAODExtension *ext = (AliAODExtension*)fExtensions->FindObject(fname);
484 ext = new AliAODExtension(fname, title);
485 fExtensions->Add(ext);
490 //______________________________________________________________________________
491 void AliAODHandler::SetOutputFileName(const char* fname)
497 //______________________________________________________________________________
498 const char *AliAODHandler::GetOutputFileName()
501 return fFileName.Data();
504 //______________________________________________________________________________
505 void AliAODHandler::SetMCHeaderInfo(AliAODMCHeader *mcHeader,AliGenEventHeader *genHeader){
508 // Utility function to cover different cases for the AliGenEventHeader
509 // Needed since different ProcessType and ImpactParamter are not
510 // in the base class...
511 // We don't encode process types for event cocktails yet
512 // could be done e.g. by adding offsets depnding on the generator
514 mcHeader->AddGeneratorName(genHeader->GetName());
515 if(!genHeader)return;
516 AliGenPythiaEventHeader *pythiaGenHeader = dynamic_cast<AliGenPythiaEventHeader*>(genHeader);
517 if (pythiaGenHeader) {
518 mcHeader->SetEventType(pythiaGenHeader->ProcessType());
522 AliGenDPMjetEventHeader* dpmJetGenHeader = dynamic_cast<AliGenDPMjetEventHeader*>(genHeader);
524 if (dpmJetGenHeader){
525 mcHeader->SetEventType(dpmJetGenHeader->ProcessType());
529 AliGenHijingEventHeader* hijingGenHeader = dynamic_cast<AliGenHijingEventHeader*>(genHeader);
531 mcHeader->SetImpactParameter(hijingGenHeader->ImpactParameter());
535 AliWarning(Form("MC Eventheader not known: %s",genHeader->GetName()));
539 ClassImp(AliAODExtension)
541 //-------------------------------------------------------------------------
542 // Support class for AOD extensions. This is created by the user analysis
543 // that requires a separate file for some AOD branches. The name of the
544 // AliAODExtension object is the file name where the AOD branches will be
546 //-------------------------------------------------------------------------
548 //______________________________________________________________________________
549 AliAODExtension::~AliAODExtension()
554 // is already handled in TerminateIO
561 //______________________________________________________________________________
562 void AliAODExtension::AddBranch(const char* cname, void* addobj)
564 // Add a new branch to the aod
567 gROOT->ProcessLine(Form("TString s_tmp; AliAnalysisManager::GetAnalysisManager()->GetAnalysisTypeString(s_tmp); sprintf((char*)0x%lx, \"%%s\", s_tmp.Data());", type));
570 TDirectory *owd = gDirectory;
574 char** apointer = (char**) addobj;
575 TObject* obj = (TObject*) *apointer;
577 fAODEvent->AddObject(obj);
579 const Int_t kSplitlevel = 99; // default value in TTree::Branch()
580 const Int_t kBufsize = 32000; // default value in TTree::Branch()
582 if (!fTreeE->FindBranch(obj->GetName())) {
583 // Do the same as if we book via
584 // TTree::Branch(TCollection*)
586 fTreeE->Bronch(obj->GetName(), cname, fAODEvent->GetList()->GetObjectRef(obj),
587 kBufsize, kSplitlevel - 1);
588 // fTreeA->Branch(obj->GetName(), cname, addobj);
593 //______________________________________________________________________________
594 Bool_t AliAODExtension::Init(Option_t *option)
597 if(!fAODEvent) fAODEvent = new AliAODEvent();
598 TDirectory *owd = gDirectory;
601 if (opt.Contains("proof")) {
603 // Merging via files. Need to access analysis manager via interpreter.
604 gROOT->ProcessLine(Form("AliAnalysisManager::GetAnalysisManager()->OpenProofFile(\"%s\", \"RECREATE\");", fName.Data()));
607 fFileE = new TFile(GetName(), "RECREATE");
609 fTreeE = new TTree("aodTree", "AliAOD tree");
610 fTreeE->Branch(fAODEvent->GetList());
615 //______________________________________________________________________________
616 Bool_t AliAODExtension::TerminateIO()