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