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 "AliAODExtension.h"
35 #include "AliAODTracklets.h"
37 #include "AliAODMCParticle.h"
38 #include "AliAODMCHeader.h"
39 #include "AliMCEventHandler.h"
40 #include "AliMCEvent.h"
41 #include "AliGenEventHeader.h"
42 #include "AliGenHijingEventHeader.h"
43 #include "AliGenDPMjetEventHeader.h"
44 #include "AliGenPythiaEventHeader.h"
45 #include "AliGenCocktailEventHeader.h"
46 #include "AliCodeTimer.h"
47 #include "AliAODBranchReplicator.h"
48 #include "Riostream.h"
50 ClassImp(AliAODHandler)
52 //______________________________________________________________________________
53 AliAODHandler::AliAODHandler() :
58 fFillExtension(kTRUE),
59 fNeedsHeaderReplication(kFALSE),
60 fNeedsTracksBranchReplication(kFALSE),
61 fNeedsVerticesBranchReplication(kFALSE),
62 fNeedsV0sBranchReplication(kFALSE),
63 fNeedsCascadesBranchReplication(kFALSE),
64 fNeedsTrackletsBranchReplication(kFALSE),
65 fNeedsPMDClustersBranchReplication(kFALSE),
66 fNeedsJetsBranchReplication(kFALSE),
67 fNeedsFMDClustersBranchReplication(kFALSE),
68 fNeedsCaloClustersBranchReplication(kFALSE),
69 fNeedsCaloTriggerBranchReplication(kFALSE),
70 fNeedsMCParticlesBranchReplication(kFALSE),
71 fNeedsDimuonsBranchReplication(kFALSE),
72 fAODIsReplicated(kFALSE),
81 // default constructor
84 //______________________________________________________________________________
85 AliAODHandler::AliAODHandler(const char* name, const char* title):
86 AliVEventHandler(name, title),
90 fFillExtension(kTRUE),
91 fNeedsHeaderReplication(kFALSE),
92 fNeedsTracksBranchReplication(kFALSE),
93 fNeedsVerticesBranchReplication(kFALSE),
94 fNeedsV0sBranchReplication(kFALSE),
95 fNeedsCascadesBranchReplication(kFALSE),
96 fNeedsTrackletsBranchReplication(kFALSE),
97 fNeedsPMDClustersBranchReplication(kFALSE),
98 fNeedsJetsBranchReplication(kFALSE),
99 fNeedsFMDClustersBranchReplication(kFALSE),
100 fNeedsCaloClustersBranchReplication(kFALSE),
101 fNeedsCaloTriggerBranchReplication(kFALSE),
102 fNeedsMCParticlesBranchReplication(kFALSE),
103 fNeedsDimuonsBranchReplication(kFALSE),
104 fAODIsReplicated(kFALSE),
113 // Normal constructor.
116 //______________________________________________________________________________
117 AliAODHandler::~AliAODHandler()
125 // is already handled in TerminateIO
135 //______________________________________________________________________________
136 Bool_t AliAODHandler::Init(Option_t* opt)
140 // Create the AODevent object
142 Bool_t createStdAOD = fIsStandard || fFillAOD;
143 if(!fAODEvent && createStdAOD){
144 fAODEvent = new AliAODEvent();
146 fAODEvent->CreateStdContent();
149 // File opening according to execution mode
153 TDirectory *owd = gDirectory;
154 if (option.Contains("proof")) {
156 // Merging via files. Need to access analysis manager via interpreter.
157 gROOT->ProcessLine(Form("AliAnalysisDataContainer *c_common_out = AliAnalysisManager::GetAnalysisManager()->GetCommonOutputContainer();"));
158 gROOT->ProcessLine(Form("AliAnalysisManager::GetAnalysisManager()->OpenProofFile(c_common_out, \"RECREATE\");"));
162 fFileA = new TFile(fFileName.Data(), "RECREATE");
168 TIter next(fExtensions);
169 AliAODExtension *ext;
170 while ((ext=(AliAODExtension*)next())) ext->Init(option);
173 TIter nextf(fFilters);
174 AliAODExtension *filteredAOD;
175 while ((filteredAOD=(AliAODExtension*)nextf())) {
176 filteredAOD->SetEvent(fAODEvent);
177 filteredAOD->Init(option);
184 //______________________________________________________________________________
185 void AliAODHandler::Print(Option_t* opt) const
187 // Print info about this object
189 cout << opt << Form("IsStandard %d filename=%s",fIsStandard,fFileName.Data()) << endl;
193 cout << opt << fExtensions->GetEntries() << " extensions :" << endl;
194 PrintExtensions(*fExtensions);
198 cout << opt << fFilters->GetEntries() << " filters :" << endl;
199 PrintExtensions(*fFilters);
203 //______________________________________________________________________________
204 void AliAODHandler::PrintExtensions(const TObjArray& array) const
206 // Show the list of aod extensions
208 AliAODExtension* ext(0x0);
209 while ( ( ext = static_cast<AliAODExtension*>(next()) ) )
215 //______________________________________________________________________________
216 void AliAODHandler::StoreMCParticles(){
219 // Remap the labels from ESD stack and store
220 // the AODMCParticles, makes only sense if we have
221 // the mcparticles branch
222 // has to be done here since we cannot know in advance
223 // which particles are needed (e.g. by the tracks etc.)
225 // Particles have been selected by AliMCEventhanlder->SelectParticle()
226 // To use the MCEventhandler here we need to set it from the outside
227 // can vanish when Handler go to the ANALYSISalice library
229 // The Branch booking for mcParticles and mcHeader has to happen
230 // in an external task for now since the AODHandler does not have access
231 // the AnalysisManager. For the same reason the pointer t o the MCEventH
232 // has to passed to the AOD Handler by this task
233 // (doing this in the steering macro would not work on PROOF)
235 if (!fAODEvent) return;
236 TClonesArray *mcarray = (TClonesArray*)fAODEvent->FindListObject(AliAODMCParticle::StdBranchName());
239 AliAODMCHeader *mcHeader = (AliAODMCHeader*)fAODEvent->FindListObject(AliAODMCHeader::StdBranchName());
242 // Get the MC Infos.. Handler needs to be set before
243 // while adding the branch
244 // This needs to be done, not to depend on the AnalysisManager
246 if(!fMCEventH)return;
247 if(!fMCEventH->MCEvent())return;
248 AliStack *pStack = fMCEventH->MCEvent()->Stack();
251 fMCEventH->CreateLabelMap();
254 // Get the Event Header
257 AliHeader* header = fMCEventH->MCEvent()->Header();
259 AliGenEventHeader* genHeader = 0;
260 if (header) genHeader = header->GenEventHeader();
263 genHeader->PrimaryVertex(vtxMC);
264 mcHeader->SetVertex(vtxMC[0],vtxMC[1],vtxMC[2]);
266 // we search the MCEventHeaders first
267 // Two cases, cocktail or not...
268 AliGenCocktailEventHeader* genCocktailHeader = dynamic_cast<AliGenCocktailEventHeader*>(genHeader);
269 if(genCocktailHeader){
270 // we have a coktail header
271 mcHeader->AddGeneratorName(genHeader->GetName());
272 // Loop from the back so that the first one sets the process type
273 TList* headerList = genCocktailHeader->GetHeaders();
274 for(int i = headerList->GetEntries()-1;i>=0;--i){
275 AliGenEventHeader *headerEntry = dynamic_cast<AliGenEventHeader*>(headerList->At(i));
277 AliFatal("AliGenEventHeader entry not found in the header list");
279 SetMCHeaderInfo(mcHeader,headerEntry);
284 // No Cocktail just take the first one
285 SetMCHeaderInfo(mcHeader,genHeader);
293 // Store the AliAODParticlesMC
294 AliMCEvent* mcEvent = fMCEventH->MCEvent();
296 Int_t np = mcEvent->GetNumberOfTracks();
297 Int_t nprim = mcEvent->GetNumberOfPrimaries();
301 TClonesArray& l = *mcarray;
303 for(int i = 0; i < np; ++i){
304 if(fMCEventH->IsParticleSelected(i)){
306 AliMCParticle* mcpart = (AliMCParticle*) mcEvent->GetTrack(i);
307 if(i<nprim)flag |= AliAODMCParticle::kPrimary;
309 if(mcEvent->IsPhysicalPrimary(i))flag |= AliAODMCParticle::kPhysicalPrim;
311 if(fMCEventH->GetNewLabel(i)!=j){
312 AliError(Form("MISMATCH New label %d j: %d",fMCEventH->GetNewLabel(i),j));
315 AliAODMCParticle mcpartTmp(mcpart,i,flag);
317 mcpartTmp.SetStatus(mcpart->Particle()->GetStatusCode());
319 Int_t d0 = mcpartTmp.GetDaughter(0);
320 Int_t d1 = mcpartTmp.GetDaughter(1);
321 Int_t m = mcpartTmp.GetMother();
323 // other than for the track labels, negative values mean
324 // no daughter/mother so preserve it
327 // no first daughter -> no second daughter
328 // nothing to be done
329 // second condition not needed just for sanity check at the end
330 mcpartTmp.SetDaughter(0,d0);
331 mcpartTmp.SetDaughter(1,d1);
332 } else if(d1 < 0 && d0 >= 0) {
334 // second condition not needed just for sanity check at the end
335 if(fMCEventH->IsParticleSelected(d0)){
336 mcpartTmp.SetDaughter(0,fMCEventH->GetNewLabel(d0));
338 mcpartTmp.SetDaughter(0,-1);
340 mcpartTmp.SetDaughter(1,d1);
342 else if (d0 > 0 && d1 > 0 ){
343 // we have two or more daughters loop on the stack to see if they are
347 for(int id = d0; id<=d1;++id){
348 if(fMCEventH->IsParticleSelected(id)){
351 d0Tmp = fMCEventH->GetNewLabel(id);
352 d1Tmp = d0Tmp; // this is to have the same schema as on the stack i.e. with one daugther d0 and d1 are the same
354 else d1Tmp = fMCEventH->GetNewLabel(id);
357 mcpartTmp.SetDaughter(0,d0Tmp);
358 mcpartTmp.SetDaughter(1,d1Tmp);
360 AliError(Form("Unxpected indices %d %d",d0,d1));
364 mcpartTmp.SetMother(m);
366 if(fMCEventH->IsParticleSelected(m))mcpartTmp.SetMother(fMCEventH->GetNewLabel(m));
367 else AliError(Form("PROBLEM Mother not selected %d \n", m));
370 new (l[j++]) AliAODMCParticle(mcpartTmp);
374 AliInfo(Form("AliAODHandler::StoreMCParticles: Selected %d (Primaries %d / total %d) after validation \n",
377 // Set the labels in the AOD output...
381 TClonesArray* tracks = fAODEvent->GetTracks();
383 for(int it = 0; it < fAODEvent->GetNTracks();++it){
384 AliAODTrack *track = fAODEvent->GetTrack(it);
387 Int_t label = track->GetLabel();
388 if(label<0){ // preserve the sign for later usage
393 if (label >= AliMCEvent::BgLabelOffset()) label = mcEvent->BgLabelToIndex(label);
394 if(label > np || track->GetLabel() == 0){
395 AliWarning(Form("Wrong ESD track label %5d (%5d)",track->GetLabel(), label));
397 if(fMCEventH->GetNewLabel(label) == 0){
398 AliWarning(Form("New label not found for %5d (%5d)",track->GetLabel(), label));
400 track->SetLabel(sign*fMCEventH->GetNewLabel(label));
405 TClonesArray *clusters = fAODEvent->GetCaloClusters();
407 for (Int_t iClust = 0;iClust < fAODEvent->GetNumberOfCaloClusters(); ++iClust) {
408 AliAODCaloCluster * cluster = fAODEvent->GetCaloCluster(iClust);
409 UInt_t nLabel = cluster->GetNLabels();
410 // Ugly but do not want to fragment memory by creating
411 // new Int_t (nLabel)
412 Int_t* labels = const_cast<Int_t*>(cluster->GetLabels());
414 for(UInt_t i = 0;i < nLabel;++i){
415 labels[i] = fMCEventH->GetNewLabel(cluster->GetLabelAt(i));
418 // cluster->SetLabels(labels,nLabel);
423 AliAODTracklets *tracklets = fAODEvent->GetTracklets();
425 for(int it = 0;it < tracklets->GetNumberOfTracklets();++it){
426 int label0 = tracklets->GetLabel(it,0);
427 int label1 = tracklets->GetLabel(it,1);
428 if(label0>=0)label0 = fMCEventH->GetNewLabel(label0);
429 if(label1>=0)label1 = fMCEventH->GetNewLabel(label1);
430 tracklets->SetLabel(it,0,label0);
431 tracklets->SetLabel(it,1,label1);
437 //______________________________________________________________________________
438 Bool_t AliAODHandler::FinishEvent()
440 // Fill data structures
441 if(fFillAOD && fFillAODRun && fAODEvent){
442 fAODEvent->MakeEntriesReferencable();
447 if ((fFillAOD && fFillAODRun) || fFillExtension) {
448 if (fExtensions && fFillExtension) {
449 // fFillExtension can be set by the ESD filter or by a delta filter in case of AOD inputs
450 TIter next(fExtensions);
451 AliAODExtension *ext;
452 while ((ext=(AliAODExtension*)next())) ext->FinishEvent();
454 if (fFilters && fFillAOD && fFillAODRun) {
455 TIter nextf(fFilters);
456 AliAODExtension *ext;
457 while ((ext=(AliAODExtension*)nextf())) {
463 if (fIsStandard && fAODEvent)
465 fAODEvent->ResetStd();
470 TClonesArray *mcarray = static_cast<TClonesArray*>(fAODEvent->FindListObject(AliAODMCParticle::StdBranchName()));
471 if(mcarray) mcarray->Delete();
473 AliAODMCHeader *mcHeader = static_cast<AliAODMCHeader*>(fAODEvent->FindListObject(AliAODMCHeader::StdBranchName()));
474 if(mcHeader) mcHeader->Reset();
477 // Reset AOD replication flag
478 fAODIsReplicated = kFALSE;
482 //______________________________________________________________________________
483 Bool_t AliAODHandler::Terminate()
486 AddAODtoTreeUserInfo();
488 TIter nextF(fFilters);
489 AliAODExtension *ext;
490 while ((ext=static_cast<AliAODExtension*>(nextF())))
492 ext->AddAODtoTreeUserInfo();
495 TIter nextE(fExtensions);
496 while ((ext=static_cast<AliAODExtension*>(nextE())))
498 ext->AddAODtoTreeUserInfo();
504 //______________________________________________________________________________
505 Bool_t AliAODHandler::TerminateIO()
513 // When closing the file, the tree is also deleted.
517 TIter nextF(fFilters);
518 AliAODExtension *ext;
519 while ((ext=static_cast<AliAODExtension*>(nextF())))
524 TIter nextE(fExtensions);
525 while ((ext=static_cast<AliAODExtension*>(nextE())))
533 //______________________________________________________________________________
534 void AliAODHandler::CreateTree(Int_t flag)
536 // Creates the AOD Tree
537 fTreeA = new TTree("aodTree", "AliAOD tree");
538 fTreeA->Branch(fAODEvent->GetList());
539 if (flag == 0) fTreeA->SetDirectory(0);
542 //______________________________________________________________________________
543 void AliAODHandler::FillTree()
550 //______________________________________________________________________________
551 void AliAODHandler::AddAODtoTreeUserInfo()
553 // Add aod event to tree user info
554 if (fTreeA) fTreeA->GetUserInfo()->Add(fAODEvent);
555 // Now the tree owns our fAODEvent...
559 //______________________________________________________________________________
560 void AliAODHandler::AddBranch(const char* cname, void* addobj, const char* filename)
562 // Add a new branch to the aod. Added optional filename parameter if the
563 // branch should be written to a separate file.
565 if (strlen(filename))
567 AliAODExtension *ext = AddExtension(filename);
568 ext->AddBranch(cname, addobj);
572 // Add branch to all filters
573 // Add branch to all filters
575 TIter next(fFilters);
576 AliAODExtension *ext;
577 while ((ext=(AliAODExtension*)next())) ext->AddBranch(cname, addobj);
580 TDirectory *owd = gDirectory;
586 char** apointer = (char**) addobj;
587 TObject* obj = (TObject*) *apointer;
589 fAODEvent->AddObject(obj);
591 const Int_t kSplitlevel = 99; // default value in TTree::Branch()
592 const Int_t kBufsize = 32000; // default value in TTree::Branch()
594 if (!fTreeA->FindBranch(obj->GetName()))
596 // Do the same as if we book via
597 // TTree::Branch(TCollection*)
599 fTreeA->Bronch(obj->GetName(), cname, fAODEvent->GetList()->GetObjectRef(obj),
600 kBufsize, kSplitlevel - 1);
605 //______________________________________________________________________________
606 AliAODExtension *AliAODHandler::AddExtension(const char *filename, const char *title)
608 // Add an AOD extension with some branches in a different file.
610 TString fname(filename);
611 if (!fname.EndsWith(".root")) fname += ".root";
613 fExtensions = new TObjArray();
614 fExtensions->SetOwner();
616 AliAODExtension *ext = (AliAODExtension*)fExtensions->FindObject(fname);
618 ext = new AliAODExtension(fname, title);
619 fExtensions->Add(ext);
624 //______________________________________________________________________________
625 AliAODExtension *AliAODHandler::GetExtension(const char *filename) const
627 // Getter for AOD extensions via file name.
628 if (!fExtensions) return NULL;
629 return (AliAODExtension*)fExtensions->FindObject(filename);
632 //______________________________________________________________________________
633 AliAODExtension *AliAODHandler::AddFilteredAOD(const char *filename, const char *filtername)
635 // Add an AOD extension that can write only AOD events that pass a user filter.
637 fFilters = new TObjArray();
638 fFilters->SetOwner();
640 AliAODExtension *filter = (AliAODExtension*)fFilters->FindObject(filename);
642 filter = new AliAODExtension(filename, filtername, kTRUE);
643 fFilters->Add(filter);
648 //______________________________________________________________________________
649 AliAODExtension *AliAODHandler::GetFilteredAOD(const char *filename) const
651 // Getter for AOD filters via file name.
652 if (!fFilters) return NULL;
653 return (AliAODExtension*)fFilters->FindObject(filename);
656 //______________________________________________________________________________
657 void AliAODHandler::SetOutputFileName(const char* fname)
663 //______________________________________________________________________________
664 const char *AliAODHandler::GetOutputFileName()
667 return fFileName.Data();
670 //______________________________________________________________________________
671 const char *AliAODHandler::GetExtraOutputs() const
673 // Get extra outputs as a string separated by commas.
674 static TString eoutputs;
678 TIter next1(fExtensions);
679 while ((obj=next1())) {
680 if (!eoutputs.IsNull()) eoutputs += ",";
681 eoutputs += obj->GetName();
685 TIter next2(fFilters);
686 while ((obj=next2())) {
687 if (!eoutputs.IsNull()) eoutputs += ",";
688 eoutputs += obj->GetName();
691 return eoutputs.Data();
694 //______________________________________________________________________________
695 Bool_t AliAODHandler::HasExtensions() const
697 // Whether or not we manage extensions
699 if ( fExtensions && fExtensions->GetEntries()>0 ) return kTRUE;
704 //______________________________________________________________________________
705 void AliAODHandler::SetMCHeaderInfo(AliAODMCHeader *mcHeader,AliGenEventHeader *genHeader){
708 // Utility function to cover different cases for the AliGenEventHeader
709 // Needed since different ProcessType and ImpactParamter are not
710 // in the base class...
711 // We don't encode process types for event cocktails yet
712 // could be done e.g. by adding offsets depnding on the generator
714 if(!genHeader)return;
715 mcHeader->AddGeneratorName(genHeader->GetName());
716 AliGenPythiaEventHeader *pythiaGenHeader = dynamic_cast<AliGenPythiaEventHeader*>(genHeader);
717 if (pythiaGenHeader) {
718 mcHeader->SetEventType(pythiaGenHeader->ProcessType());
719 mcHeader->SetPtHard(pythiaGenHeader->GetPtHard());
723 AliGenDPMjetEventHeader* dpmJetGenHeader = dynamic_cast<AliGenDPMjetEventHeader*>(genHeader);
725 if (dpmJetGenHeader){
726 mcHeader->SetEventType(dpmJetGenHeader->ProcessType());
730 AliGenHijingEventHeader* hijingGenHeader = dynamic_cast<AliGenHijingEventHeader*>(genHeader);
732 mcHeader->SetImpactParameter(hijingGenHeader->ImpactParameter());
736 AliWarning(Form("MC Eventheader not known: %s",genHeader->GetName()));