]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGGA/GammaConv/AliAnalysisTaskdPhi.cxx
0c3f696483597f7825e9f4796f53e535598bc1b3
[u/mrichter/AliRoot.git] / PWGGA / GammaConv / AliAnalysisTaskdPhi.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                                                                                        *
4  * Authors: Svein Lindal                                                  *
5  * Version 1.0                                                                    *
6  *                                                                                                                                        *
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 ////////////////////////////////////////////////
18 //--------------------------------------------- 
19 // Class doing conversion gamma dPhi correlations
20 // Gamma Conversion analysis
21 //---------------------------------------------
22 ////////////////////////////////////////////////
23
24 #include "AliAnalysisTaskdPhi.h"
25
26 #include <TH2I.h>
27 #include <TList.h>
28 #include <TChain.h>
29
30 #include <AliAnalysisManager.h>
31 #include <AliInputEventHandler.h>
32 #include <AliESDInputHandler.h>
33 #include <AliAODInputHandler.h>
34
35 #include "AliConversionTrackCuts.h"
36 #include "AliConversionCuts.h"
37 #include "AliConversionMesonCuts.h"
38 #include "AliAODConversionPhoton.h"
39 #include "AliAODConversionMother.h"
40 #include "AliAnaConvCorrPhoton.h"
41 #include "AliAnaConvCorrPion.h"
42 #include "AliAnaConvIsolation.h"
43 #include "AliV0ReaderV1.h"
44 // Author Svein Lindal <slindal@fys.uio.no>
45 using namespace std;
46
47 ClassImp(AliAnalysisTaskdPhi)
48
49
50 //________________________________________________________________________
51 AliAnalysisTaskdPhi::AliAnalysisTaskdPhi(const char *name) : AliAnalysisTaskSE(name),
52   fHistograms(NULL),
53   fV0Reader(NULL),
54   fSaveReaderHists(kFALSE),
55   fV0Filter(NULL),
56   fPhotonFilter(NULL),
57   fMesonFilter(NULL),
58   fTrackCuts(NULL),
59   hMEvents(NULL),
60   hTrackCent(NULL),
61   fPhotonCorr(NULL),
62   fPionCorr(NULL), 
63   fDeltaAODBranchName("AliAODGammaConversion_gamma"), 
64   fAxistPt(),
65   fAxiscPt(),
66   fAxisEta(),
67   fAxisPhi(),
68   fAxisCent(),
69   fAxisZ(), 
70   fAxisPiM(), 
71   fDoPhoton(kFALSE)
72 {
73   //constructor
74   fAxistPt.SetNameTitle("tPtAxis", "tPt");
75   fAxistPt.Set(20, 0, 100);
76
77   fAxiscPt.SetNameTitle("cPtAxis", "cPt");
78   fAxiscPt.Set(20, 0, 100);
79
80   fAxisEta.SetNameTitle("EtaAxis", "Eta");
81   fAxisEta.Set(180, -0.9, 0.9);
82
83   fAxisPhi.SetNameTitle("PhiAxis", "Phi");
84   fAxisPhi.Set(128, 0, TMath::TwoPi());
85
86   fAxisZ.SetNameTitle("ZAxis", "Z");
87   fAxisZ.Set(4, -10, 10);
88
89   fAxisCent.SetNameTitle("CentAxis", "Cent");
90
91   Double_t centbins[5] = {0, 10, 30, 60, 100.1};
92   fAxisCent.Set(4, centbins);
93
94   Double_t mbins[7] = {0.1, 0.11, 0.12, 0.15, 0.16, 0.18, 0.2};
95   fAxisPiM.SetNameTitle("InvMassPi0", "Invariant mass");
96   fAxisPiM.Set(6, mbins);
97
98
99   DefineInput(0, TChain::Class());
100   DefineOutput(1, TList::Class());
101 }
102
103
104
105 //________________________________________________________________________
106 AliAnalysisTaskdPhi::~AliAnalysisTaskdPhi(){
107   //destructor
108   
109   if(fV0Filter)
110         delete fV0Filter;
111   fV0Filter = NULL;
112
113   if(fMesonFilter)
114         delete fMesonFilter;
115   fMesonFilter = NULL;
116
117   if(fPhotonFilter)
118         delete fPhotonFilter;
119   fPhotonFilter = NULL;
120
121   if(fHistograms)
122         delete fHistograms;
123   fHistograms = NULL;
124
125   if(fTrackCuts)
126     delete fTrackCuts;
127   fTrackCuts = NULL;
128
129 }
130
131 ///________________________________________________________________________
132 void AliAnalysisTaskdPhi::SetUpCorrObjects() {
133   ///Creat corr obj
134   //  fIsoAna = new AliAnaConvIsolation();
135   
136   AliDebug(AliLog::kDebug + 5, "Set Up corr objects");
137
138   
139
140
141   if(fDoPhoton) {
142     fPhotonCorr = new AliAnaConvCorrPhoton("PhotonCorr","photon %s");
143     fPhotonCorr->GetAxisCent().Set(fAxisCent.GetNbins(), fAxisCent.GetXbins()->GetArray());
144     fPhotonCorr->GetAxisZ().Set(fAxisZ.GetNbins(), fAxisZ.GetXbins()->GetArray());
145     fPhotonCorr->GetAxistPt().Set(fAxistPt.GetNbins(), fAxistPt.GetXbins()->GetArray());
146     fPhotonCorr->GetAxiscPt().Set(fAxiscPt.GetNbins(), fAxiscPt.GetXbins()->GetArray());
147     fPhotonCorr->CreateHistograms();
148     fHistograms->Add(fPhotonCorr->GetHistograms());
149   }
150   
151   fPionCorr = new AliAnaConvCorrPion("PionCorr", "pion");
152   fPionCorr->GetAxisCent().Set(fAxisCent.GetNbins(), fAxisCent.GetXbins()->GetArray());
153   fPionCorr->GetAxisZ().Set(fAxisZ.GetNbins(), fAxisZ.GetXbins()->GetArray());
154   fPionCorr->GetAxistPt().Set(fAxistPt.GetNbins(), fAxistPt.GetXbins()->GetArray());
155   fPionCorr->GetAxiscPt().Set(fAxiscPt.GetNbins(), fAxiscPt.GetXbins()->GetArray());
156   fPionCorr->GetAxisM().Set(fAxisPiM.GetNbins(), fAxisPiM.GetXbins()->GetArray());
157   fPionCorr->CreateHistograms();
158   fHistograms->Add(fPionCorr->GetHistograms());
159 }
160
161     
162     
163         
164           
165
166
167
168 //________________________________________________________________________
169 void AliAnalysisTaskdPhi::UserCreateOutputObjects() {
170   // Create histograms
171   
172   fHistograms = new TList();
173   fHistograms->SetName("dPhi_histograms");
174   fHistograms->SetOwner(kTRUE);
175
176
177   if(!fV0Reader){
178     fV0Reader=(AliV0ReaderV1*)AliAnalysisManager::GetAnalysisManager()->GetTask("V0ReaderV1");
179   }
180   
181   if(!fV0Reader){
182     printf("Error: No V0 Reader");
183   } // GetV0Reader
184
185
186   if(fSaveReaderHists) {
187     AliConversionCuts * v0cuts = fV0Reader->GetConversionCuts();
188     if(v0cuts) {
189       TList * histograms = v0cuts->GetCutHistograms();
190       if(!histograms) {
191         AliWarning("initializing v0 reader hists");
192         v0cuts->InitCutHistograms("V0Reader", kTRUE);
193       }
194       histograms = v0cuts->GetCutHistograms();
195       if(histograms) {
196         fHistograms->Add(histograms);
197       }
198     }
199   }
200
201   if(fV0Filter) {
202     fV0Filter->InitCutHistograms("V0Filter", kFALSE);
203     fHistograms->Add(fV0Filter->GetCutHistograms());
204   }
205   if(fMesonFilter) {
206     fMesonFilter->InitCutHistograms("PionFilter", kFALSE);
207     fHistograms->Add(fMesonFilter->GetCutHistograms());
208   }
209   if(fPhotonFilter) {
210     fPhotonFilter->InitCutHistograms("PhotonFilter", kFALSE);
211     fHistograms->Add(fPhotonFilter->GetCutHistograms());
212   }
213   
214
215   AliConversionTrackCuts * tc = dynamic_cast<AliConversionTrackCuts*>(fTrackCuts);
216   if(tc) fHistograms->Add(tc->CreateHistograms());
217
218   SetUpCorrObjects();
219
220
221   ///Set up ME histograms
222   TList * MEHistograms = new TList();
223   MEHistograms->SetName("MEHistograms");
224   MEHistograms->SetOwner(kTRUE);
225   fHistograms->Add(MEHistograms);
226
227   // hMETracks = new TObjArray();
228   // hMETracks->SetName("TrackArray");
229   // hMETracks->SetOwner(kTRUE);
230   // hMEPhotons = new TObjArray();
231   // hMEPhotons->SetName("PhotonArray");
232   // hMEPhotons->SetOwner(kTRUE);
233   // hMEPions = new TObjArray();
234   // hMEPions->SetName("PionArray");
235   // hMEPions->SetOwner(kTRUE);
236
237   // MEHistograms->Add(hMETracks);
238   // MEHistograms->Add(hMEPions);
239   // MEHistograms->Add(hMEPhotons);
240
241   hMEvents = new TH2I("hMEvents", "Nevents vs centrality vertexz",
242                       fAxisZ.GetNbins(), fAxisZ.GetXbins()->GetArray(),
243                       fAxisCent.GetNbins(), fAxisCent.GetXbins()->GetArray());
244   MEHistograms->Add(hMEvents);
245
246   hTrackCent = new TH2I("hTrackCent", "N accepted tracks vs centrality",
247                         fAxisCent.GetNbins() > 2 ? 100 : 1, -999, 999,
248                         750, 0, 1500);
249   MEHistograms->Add(hTrackCent);
250
251   // TList axesList;
252   // axesList.AddAt(&GetAxisEta(), 0);
253   // axesList.AddAt(&GetAxisPhi(), 1);
254   // axesList.AddAt(&GetAxistPt(), 2);
255   // axesList.SetOwner(kFALSE);
256   
257   // TList piAxesList;
258   // piAxesList.AddAt(&GetAxisEta(), 0);
259   // piAxesList.AddAt(&GetAxisPhi(), 1);
260   // piAxesList.AddAt(&GetAxistPt(), 2);
261   // piAxesList.AddAt(&GetAxisPiMass(), 3);
262   // piAxesList.SetOwner(kFALSE);
263
264   // TList * outAxesList = new TList();
265   // outAxesList->Add(&fAxisCent);
266   // outAxesList->Add(&fAxisZ);
267   // fHistograms->Add(outAxesList);
268
269   // for(Int_t iz = 0; iz < fAxisZ.GetNbins(); iz++) {
270   //    TObjArray * trackArray = new TObjArray();
271   //    trackArray->SetName(Form("METracks_%d", iz));
272   //    trackArray->SetOwner(kTRUE);
273   //    TObjArray * photonArray = new TObjArray();
274   //    photonArray->SetName(Form("MEPhotons_%d", iz));
275   //    photonArray->SetOwner(kTRUE);
276   //    TObjArray * pionArray = new TObjArray();
277   //    pionArray->SetName(Form("MEPions_%d", iz));
278   //    pionArray->SetOwner(kTRUE);
279
280
281   //    hMEPions->AddAt(pionArray, iz);
282   //    hMETracks->AddAt(trackArray, iz);
283   //    hMEPhotons->AddAt(photonArray, iz);
284
285   //    for(Int_t ic = 0; ic < fAxisCent.GetNbins(); ic++) {
286
287   //      TString nameString = Form("%d_%d", iz, ic);
288   //      TString titleString = Form("%f < Z < %f ... %f cent %f", 
289   //                                                             fAxisZ.GetBinLowEdge(iz+1), fAxisZ.GetBinUpEdge(iz+1), 
290   //                                                             fAxisCent.GetBinLowEdge(ic+1), fAxisCent.GetBinUpEdge(ic+1));
291
292
293   //      THnSparseF * trackHistogram = CreateSparse(Form("tracks_%s", nameString.Data()), 
294   //                                                                                             Form("tracks %s", titleString.Data()), &axesList );
295   //      trackArray->AddAt(trackHistogram, ic);
296
297   //      THnSparseF * photonHistogram = CreateSparse(Form("photons_%s", nameString.Data()), 
298   //                                                                                             Form("photons %s", titleString.Data()), &axesList );
299   //      photonArray->AddAt(photonHistogram, ic);
300
301   //      THnSparseF * pionHistogram = CreateSparse(Form("pions_%s", nameString.Data()), 
302   //                                                                                             Form("pions %s", titleString.Data()), &piAxesList );
303   //      pionArray->AddAt(pionHistogram, ic);
304   //    }
305   // }
306
307   PostData(1, fHistograms);
308
309 }
310
311 ///________________________________________________________________________
312 THnSparseF * AliAnalysisTaskdPhi::CreateSparse(TString nameString, TString titleString, TList * axesList) {
313   ///Create sparse
314   const Int_t dim = axesList->GetSize();
315
316   TAxis * axes[dim];
317   Int_t bins[dim];
318   Double_t min[dim];
319   Double_t max[dim];
320
321   for(Int_t i = 0; i<dim; i++) {
322         TAxis * axis = dynamic_cast<TAxis*>(axesList->At(i));
323         if(axis) {
324           axes[i] = axis;
325         } else {
326           cout << "AliAnalysisTaskdPhi::CreateSparse: Error error, all the axes are not present in axis list" << endl;
327           return NULL;
328         }
329   }
330
331   for(Int_t i = 0; i<dim; i++) {
332         bins[i] = axes[i]->GetNbins(); 
333         min[i] = axes[i]->GetBinLowEdge(1);
334         max[i] = axes[i]->GetBinUpEdge(axes[i]->GetNbins());
335   }
336
337   THnSparseF * sparse = new THnSparseF(Form("METracks_%s", nameString.Data()), 
338                                                                                            Form("tracks %s", titleString.Data()), 
339                                                                                            dim, bins, min, max);
340   
341   for(Int_t i = 0; i<dim; i++) {
342         sparse->GetAxis(i)->SetNameTitle(axes[i]->GetName(), axes[i]->GetTitle() );
343         if(axes[i]->GetXbins()->GetSize() > 0) {
344           sparse->SetBinEdges(i, axes[i]->GetXbins()->GetArray() );
345         }
346   }
347
348   return sparse;
349 }
350
351 //________________________________________________________________________
352 void AliAnalysisTaskdPhi::UserExec(Option_t *) {
353   ///User exec. 
354
355   //if(! fV0Filter->EventIsSelected(fInputEvent)) return;
356   if(!fV0Reader){
357     AliError("Error: No V0 Reader");
358     return;
359   } // GetV0Reader
360
361   if(!fV0Reader->IsEventSelected()) {
362         return;
363   }
364    AliDebug(5, "Processing event");
365  
366   AliAnalysisManager *man=AliAnalysisManager::GetAnalysisManager();
367   Bool_t isAOD=man->GetInputEventHandler()->IsA()==AliAODInputHandler::Class();
368   
369   AliInputEventHandler* inputHandler = (AliInputEventHandler*) (man->GetInputEventHandler());
370   if (!inputHandler) {
371         cout << "cout no input event handler"<<endl;
372         return;
373   }
374
375   
376   if ( fV0Filter && !fV0Filter->GetPIDResponse() ) {
377     cout << "aaaa"<<endl;
378         if ( inputHandler->GetPIDResponse() ){
379           fV0Filter->SetPIDResponse( inputHandler->GetPIDResponse() );
380         } else {
381           
382           //AOD case
383           if (isAOD){
384                 if (!fV0Filter->GetPIDResponse()){
385                   fV0Filter->InitAODpidUtil(1);
386                 }
387           }
388         }
389   }
390
391   Double_t centrality = 0.0;
392   Double_t eventPlane = 0.0;
393   Double_t vertexz = fInputEvent->GetPrimaryVertex()->GetZ();
394   if(isAOD) {
395     AliAODHeader * header = static_cast<AliAODHeader*>(fInputEvent->GetHeader());
396         centrality = header->GetCentrality();
397         eventPlane = header->GetEventplane();
398   } else {
399         centrality = static_cast<AliESDEvent*>(fInputEvent)->GetCentrality()->GetCentralityPercentile("V0M");
400         eventPlane = fInputEvent->GetEventplane()->GetEventplane("Q");
401   }
402
403
404   const Int_t centBin = GetBin(fAxisCent, centrality);
405   const Int_t vertexBin = GetBin(fAxisZ, vertexz);
406
407
408   if(DebugLevel () > 4) {
409         cout << "centrality: " << centrality <<  " " << GetBin(fAxisCent, centrality) << endl;
410         cout << "vertexz: " << vertexz <<  " " << GetBin(fAxisZ, vertexz) << endl;
411         cout << "eventPlane: " << eventPlane <<  " " << endl;
412   }
413
414
415
416   if(centBin < 0 || vertexBin < 0) {
417     //  AliError("bin out of range");
418         return;
419   }
420
421   TObjArray gammas;
422
423   //TClonesArray * aodGammas = GetConversionGammas(isAOD);
424   TClonesArray * aodGammas = fV0Reader->GetReconstructedGammas();
425   if(!aodGammas) {
426     AliError("no aod gammas found!");
427     return;
428   }
429
430   
431   if(DebugLevel() > 1) printf("Number of conversion gammas %d \n", aodGammas->GetEntriesFast());
432   for(Int_t ig = 0; ig < aodGammas->GetEntriesFast(); ig++) {
433     AliAODConversionPhoton * photon = dynamic_cast<AliAODConversionPhoton*>(aodGammas->At(ig));
434     
435     if(!photon) continue;
436     if(!fV0Filter || fV0Filter->PhotonIsSelected(static_cast<AliConversionPhotonBase*>(photon), fInputEvent)) {
437       gammas.Add(static_cast<TObject*>(photon));
438     }
439   }
440   
441   if(DebugLevel() > 4) printf("Number of accepted gammas %d \n", gammas.GetEntriesFast());
442   hMEvents->Fill(vertexz, centrality);
443   
444   ///create track array
445   TObjArray tracks;
446   const Double_t etalim[2] = { fAxisEta.GetBinLowEdge(1), fAxisEta.GetBinUpEdge(fAxisEta.GetNbins())};
447   for(Int_t iTrack = 0; iTrack < fInputEvent->GetNumberOfTracks(); iTrack++) {
448     AliVTrack * track = static_cast<AliVTrack*>(fInputEvent->GetTrack(iTrack));
449     if(track->Pt() < fAxiscPt.GetBinLowEdge(1) ) continue;
450     if(track->Eta() < etalim[0] || track->Eta() > etalim[1]) continue;
451     if(!fTrackCuts || fTrackCuts->IsSelected(track)) {
452       tracks.Add(track);
453     }
454   }
455
456   hTrackCent->Fill(centrality, tracks.GetEntriesFast());
457   Process(&gammas, &tracks, centrality, vertexz);
458
459   PostData(1, fHistograms);
460   
461 }
462
463 //________________________________________________________________________
464 void AliAnalysisTaskdPhi::Process(TObjArray * gammas, TObjArray * tracks, const Float_t cent, const Float_t vtxz) {
465   ///Process stuff
466
467   const Double_t etalim[2] = { fAxisEta.GetBinLowEdge(1), fAxisEta.GetBinUpEdge(fAxisEta.GetNbins())};
468   if(DebugLevel() > 4) printf("Number of accepted gammas, tracks %d  %d \n", gammas->GetEntriesFast(), tracks->GetEntriesFast());
469  
470   AliAnaConvCorrBase * gCorr = fPhotonCorr; //GetCorrObject(vertexBin, centBin, fPhotonCorr);
471   AliAnaConvCorrPion * piCorr = fPionCorr; //static_cast<AliAnaConvCorrPion*>(GetCorrObject(vertexBin, centBin, fPionCorr));
472   
473   if(!piCorr) {
474     AliError("corr object missing");
475     return;
476   }
477   
478   TObjArray photons;
479   TObjArray pions;
480
481   for(Int_t i1 = 0; i1 < gammas->GetEntriesFast(); i1++) {
482     AliAODConversionPhoton * ph1 = static_cast<AliAODConversionPhoton*>(gammas->UncheckedAt(i1));
483     Int_t tIDs[4] = {ph1->GetLabel(0), ph1->GetLabel(1), -1, -1};
484         
485
486     if(gCorr) {
487       if(!fPhotonFilter || fPhotonFilter->PhotonIsSelected(static_cast<AliConversionPhotonBase*>(ph1), fInputEvent)) {
488         if(ph1->Pt() > fAxistPt.GetBinLowEdge(1)) {
489           if(ph1->Eta() > etalim[0] && ph1->Eta() < etalim[1]) {
490             gCorr->CorrelateWithTracks( static_cast<AliAODConversionParticle*>(ph1), tracks, tIDs, cent, vtxz);
491             photons.Add(ph1);
492           }
493         }
494       }
495     }
496
497     for(Int_t i2 = 0; i2 < i1; i2++) {
498       AliAODConversionPhoton * ph2 = static_cast<AliAODConversionPhoton*>(gammas->UncheckedAt(i2));
499           
500       if( ph2->GetTrackLabelPositive()==ph1->GetTrackLabelPositive() 
501           || ph2->GetTrackLabelNegative()==ph1->GetTrackLabelNegative()
502           || ph2->GetTrackLabelNegative()==ph1->GetTrackLabelPositive()
503           || ph2->GetTrackLabelPositive()==ph1->GetTrackLabelNegative()) {
504         continue;
505       }
506
507       AliAODConversionMother * pion = new AliAODConversionMother(ph1, ph2);
508       pion->SetLabels(i1, i2);
509       pion->CalculateDistanceOfClossetApproachToPrimVtx(fInputEvent->GetPrimaryVertex());
510       
511       if(!fMesonFilter || fMesonFilter->MesonIsSelected(pion, kTRUE) ) {
512         tIDs[2] = ph2->GetLabel(0);
513         tIDs[3] = ph2->GetLabel(1);
514         piCorr->FillTriggerCounters(pion);
515         AliDebug(AliLog::kDebug + 5, "We have a pion");
516         if(pion->Pt() > fAxistPt.GetBinLowEdge(1) && 
517            pion->M() > fAxisPiM.GetBinLowEdge(1) && 
518            pion->M() < fAxisPiM.GetBinUpEdge(fAxisPiM.GetNbins()) &&
519            pion->Eta() > etalim[0] && pion->Eta() < etalim[1]) {
520           piCorr->CorrelateWithTracks(pion, tracks, tIDs, cent, vtxz);
521           pions.Add(static_cast<TObject*>(pion));
522         }
523       }
524     }
525   }
526
527   piCorr->FillCounters(&pions, tracks, cent, vtxz);
528   if(gCorr) {
529     gCorr->FillCounters(&photons, tracks, cent, vtxz);
530   }
531 }
532
533 //________________________________________________________________________
534 void AliAnalysisTaskdPhi::Terminate(Option_t *) {
535  
536   // Draw result to the screen
537   // Called once at the end of the query
538 }
539
540 //________________________________________________________________________
541 TClonesArray * AliAnalysisTaskdPhi::GetConversionGammas(Bool_t isAOD) {
542   
543
544   if(isAOD) {
545
546         TClonesArray * gammas = dynamic_cast<TClonesArray*>(fInputEvent->FindListObject(fDeltaAODBranchName.Data()));
547         if(gammas) {
548           return gammas;
549         }
550
551         FindDeltaAODBranchName(fInputEvent);
552         gammas = dynamic_cast<TClonesArray*>(fInputEvent->FindListObject(fDeltaAODBranchName.Data()));
553         return gammas;
554
555   } else {
556     TClonesArray * gammas = dynamic_cast<TClonesArray*>(GetInputData(1));
557     return gammas;
558   }
559
560 }
561
562 //________________________________________________________________________
563 void AliAnalysisTaskdPhi::FindDeltaAODBranchName(AliVEvent * event){
564   ///Find aod branch
565   TList *list=event->GetList();
566   for(Int_t ii=0;ii<list->GetEntries();ii++){
567         TString name((list->At(ii))->GetName());
568         if(name.BeginsWith("GammaConv")&&name.EndsWith("gamma")){
569           fDeltaAODBranchName=name;
570           AliDebug(AliLog::kDebug + 5, Form("Set DeltaAOD BranchName to: %s",fDeltaAODBranchName.Data()));
571           return;
572         }
573   }
574 }
575   
576