]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/Correlations/DPhi/AliAnalysisTaskLongRangeCorrelations.cxx
coverity fixes + adding HM toy model in CMakelibPWGCFebye.pkg and PWGCFebyeLinkDef.h
[u/mrichter/AliRoot.git] / PWGCF / Correlations / DPhi / AliAnalysisTaskLongRangeCorrelations.cxx
1 // -*- C++ -*-
2 /**************************************************************************
3  * Copyright(c) 2012,      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 // $Id: AliAnalysisTaskLongRangeCorrelations.cxx 239 2012-12-12 17:25:09Z cmayer $
17
18 #include <numeric>
19 #include <functional>
20
21 #include <TChain.h>
22 #include <THashList.h>
23 #include <THnSparse.h>
24 #include <TH1.h>
25 #include <TH2.h>
26 #include <TH3.h>
27
28 #include "AliEventPoolManager.h"
29 #include "AliAnalysisManager.h"
30 #include "AliAODEvent.h"
31 #include "AliAODMCHeader.h"
32 #include "AliAODMCParticle.h"
33 #include "AliMCEvent.h"
34 #include "AliInputEventHandler.h"
35 #include "AliLog.h"
36
37 #include "AliAnalysisTaskLongRangeCorrelations.h"
38
39
40 ClassImp(AliAnalysisTaskLongRangeCorrelations);
41 ClassImp(LRCParticle);
42
43 AliAnalysisTaskLongRangeCorrelations::AliAnalysisTaskLongRangeCorrelations(const char *name)
44   : AliAnalysisTaskSE(name)
45   , fOutputList(NULL)
46   , fRunMixing(kFALSE)
47   , fPoolMgr(NULL)
48   , fMixingTracks(50000)
49   , fTrackFilter(128)
50   , fCentMin(0), fCentMax(20)
51   , fPtMin(0.2), fPtMax(1e10)
52   , fPhiMin(0.), fPhiMax(TMath::TwoPi())
53   , fMaxAbsVertexZ(10.)
54   , fnBinsCent( 220), fnBinsPt(400), fnBinsPhi(4),              fnBinsEta(120)
55   , fxMinCent( -5.0), fxMinPt( 0.0), fxMinPhi( 0.0),            fxMinEta(-1.5)
56   , fxMaxCent(105.0), fxMaxPt( 4.0), fxMaxPhi( TMath::TwoPi()), fxMaxEta( 1.5) {
57   DefineInput(0, TChain::Class());
58   DefineOutput(1, TList::Class());
59 }
60
61 AliAnalysisTaskLongRangeCorrelations::~AliAnalysisTaskLongRangeCorrelations() {
62   if (NULL != fOutputList) {
63     delete fOutputList;
64     fOutputList = NULL;
65   }
66 }
67
68 void AliAnalysisTaskLongRangeCorrelations::UserCreateOutputObjects() {
69   fOutputList = new THashList;
70   fOutputList->SetOwner(kTRUE);
71   fOutputList->SetName(GetOutputListName());
72
73   // Event statistics
74   const char* eventStatLabels[] = { 
75     "All Events",
76     "Physics Selection",
77     "Centrality Selection",
78     "Vertex Selection",
79     "Analyzed Events"
80   };
81   const size_t nEventStat(sizeof(eventStatLabels)/sizeof(const char*));
82   TH1* hStats(new TH1D("histEventStats", "histEventStats", nEventStat, -0.5, nEventStat-0.5));
83   for (size_t i(0); i<nEventStat; ++i)
84     hStats->GetXaxis()->SetBinLabel(1+i, eventStatLabels[i]);
85   fOutputList->Add(hStats);
86
87   // QA histograms
88   fOutputList->Add(new TH1D("histQAVertexZ", ";vertex-z (cm)",
89                             800, -40, 40));
90   fOutputList->Add(new TH3D("histQACentPt", ";charge;centrality V0M(%);p_{T} (GeV/c)",
91                             2, -0.5, 1.5, fnBinsCent, fxMinCent, fxMaxCent, fnBinsPt, fxMinPt, fxMaxPt));
92   fOutputList->Add(new TH3D("histQAPhiEta", ";charge;#phi (rad);#eta",
93                             2, -0.5, 1.5, 200, 0.0, TMath::TwoPi(), 300, -1.5, 1.5));
94
95   // N(eta) distributions with different binnings
96   fOutputList->Add(new TH2D("histNEta_300", ";#eta;N", 300, -1.5, 1.5, 1000, 0., 1000.));  // 0.01
97   fOutputList->Add(new TH2D("histNEta_120", ";#eta;N", 120, -1.5, 1.5, 1000, 0., 1000.));  // 0.025
98   fOutputList->Add(new TH2D("histNEta__30", ";#eta;N",  30, -1.5, 1.5, 1000, 0., 1000.));  // 0.1
99   fOutputList->Add(new TH2D("histNEta__15", ";#eta;N",  15, -1.5, 1.5, 1000, 0., 1000.));  // 0.2
100
101   // Moments
102   fOutputList->Add(MakeHistSparsePhiEta("histMoment1PhiEta_1"));
103   fOutputList->Add(MakeHistSparsePhiEta("histMoment1PhiEta_2"));
104   fOutputList->Add(MakeHistSparsePhiEtaPhiEta("histMoment2PhiEtaPhiEta_11"));
105   fOutputList->Add(MakeHistSparsePhiEtaPhiEta("histMoment2PhiEtaPhiEta_12"));
106   fOutputList->Add(MakeHistSparsePhiEtaPhiEta("histMoment2PhiEtaPhiEta_22"));
107   
108   // add MC Histograms
109   const Int_t N(fOutputList->GetEntries());
110   for (Int_t i(0); i<N; ++i)
111     fOutputList->Add(fOutputList->At(i)->Clone(TString("MC_")
112                                                + fOutputList->At(i)->GetName()));
113
114   if (fRunMixing)
115     SetupForMixing();
116
117   PostData(1, fOutputList);
118 }
119
120 void AliAnalysisTaskLongRangeCorrelations::UserExec(Option_t* ) {
121   AliAnalysisManager* pManager(AliAnalysisManager::GetAnalysisManager());
122   if (NULL == pManager) return;
123
124   AliInputEventHandler* pInputHandler(dynamic_cast<AliInputEventHandler*>(pManager->GetInputEventHandler()));
125   if (NULL == pInputHandler) return;
126
127   AliAODEvent* pAOD(dynamic_cast<AliAODEvent*>(InputEvent()));
128   if (NULL == pAOD) return;
129
130   AliAODHeader *pAODHeader = pAOD->GetHeader();
131   if (NULL == pAODHeader) return;
132
133   AliAODMCHeader *pAODMCHeader(dynamic_cast<AliAODMCHeader*>(pAOD->FindListObject(AliAODMCHeader::StdBranchName())));
134   const Bool_t isMC(NULL != pAODMCHeader);
135
136   TClonesArray *arrayMC(NULL);
137   if (isMC) {
138     arrayMC = dynamic_cast<TClonesArray*>(pAOD->GetList()->FindObject(AliAODMCParticle::StdBranchName()));
139     if (NULL == arrayMC) return;
140   }
141
142   // --- event cuts MC ---
143   if (isMC) {
144     Fill("MC_histEventStats", 0.); // all events
145     Fill("MC_histEventStats", 1.); // events passing physics selection
146     Fill("MC_histEventStats", 2.); // events passing centrality selection
147     Fill("MC_histQAVertexZ", pAODMCHeader->GetVtxZ());
148
149     // vertex cut in MC
150     if (TMath::Abs(pAODMCHeader->GetVtxZ()) > fMaxAbsVertexZ) return;    
151     Fill("MC_histEventStats", 3.); // events passing vertex selection
152   }
153
154   // --- event cuts data ---
155   Fill("histEventStats", 0.); // all events
156
157   if (!pInputHandler->IsEventSelected()) return;
158   Fill("histEventStats", 1.); // events passing physics selection
159
160   const Double_t centrality(pAODHeader->GetCentralityP()->GetCentralityPercentile("V0M"));
161   if (centrality < fCentMin || centrality >= fCentMax) return;
162   Fill("histEventStats", 2.); // events passing centrality selection
163
164   // vertex selection
165   const Int_t nVertex(pAOD->GetNumberOfVertices());
166   if (0 == nVertex) return;
167   const AliAODVertex* pVertex(pAOD->GetPrimaryVertex());
168   if (NULL == pVertex) return;
169   const Int_t nTracksPrimary(pVertex->GetNContributors());
170   if (nTracksPrimary < 1) return;
171
172   const Double_t zVertex(pVertex->GetZ());
173   Fill("histQAVertexZ", zVertex);
174   if (TMath::Abs(zVertex) > fMaxAbsVertexZ) return;
175
176   Fill("histEventStats", 3.); // events passing vertex selection
177
178   // event is accepted
179   TObjArray* tracksMain(GetAcceptedTracks(pAOD, centrality));
180
181   if (fRunMixing) {  
182     AliEventPool* pEventPool(fPoolMgr->GetEventPool(centrality, pAOD->GetPrimaryVertex()->GetZ()));
183     if (NULL == pEventPool)
184       AliFatal(Form("No pool found for centrality = %f, zVtx = %f", centrality, pAOD->GetPrimaryVertex()->GetZ()));
185     
186 //     pEventPool->PrintInfo();
187     if (pEventPool->IsReady()
188         || pEventPool->NTracksInPool()     > fMixingTracks/10
189         || pEventPool->GetCurrentNEvents() >= 5) {
190       Fill("histEventStats", 4.); // analyzed events
191       const Int_t nMix(pEventPool->GetCurrentNEvents());
192       for (Int_t i(0); i<nMix; ++i) {
193         TObjArray* tracksMixed(pEventPool->GetEvent(i));
194         CalculateMoments("", tracksMain, tracksMixed, 1./nMix);
195       }
196     }
197     // Update the Event pool
198     pEventPool->UpdatePool(tracksMain);
199   } else { // no mixing
200     Fill("histEventStats", 4.); // analyzed events
201     CalculateMoments("", tracksMain, tracksMain, 1.);
202     delete tracksMain;
203   }
204   
205   if (isMC) {
206     TObjArray* tracksMC(GetAcceptedTracks(arrayMC, centrality));
207     Fill("MC_histEventStats", 4.); // analyzed MC events
208     CalculateMoments("MC_", tracksMC, tracksMC, 1.);
209     delete tracksMC;
210   }
211 }
212
213 void AliAnalysisTaskLongRangeCorrelations::Terminate(Option_t* ) {
214   //
215 }
216
217 TString AliAnalysisTaskLongRangeCorrelations::GetOutputListName() const {
218   TString listName("listLRC");
219   listName += TString::Format("_%smix",         fRunMixing ? "" : "no");
220   listName += TString::Format("_trackFilt%d",   fTrackFilter);
221   listName += TString::Format("_cent%.0fT%.0f", fCentMin, fCentMax);
222   listName += TString::Format("_ptMin%.0fMeV",  1e3*fPtMin);
223   listName += TString::Format("_phi%.0fT%.0f",  TMath::RadToDeg()*fPhiMin, TMath::RadToDeg()*fPhiMax);
224   return listName;
225 }
226
227 void AliAnalysisTaskLongRangeCorrelations::SetupForMixing() {
228   const Int_t trackDepth(fMixingTracks);
229   const Int_t poolsize(1000); // Maximum number of events
230
231   Double_t centralityBins[] = { // centrality bins
232 //     0.,1.,2.,3.,4.,5.,6.,7.,8.,9.,10.,15.,20.,25.,30.,35.,40.,45.,50.,55.,60.,65.,70.,75.,80.,90.,100.
233      0.,20.,100.
234   };
235   const Int_t nCentralityBins(sizeof(centralityBins)/sizeof(Double_t) - 1);
236
237   Double_t vertexBins[] = {  // Zvtx bins
238     -10., -7., -5., -3., -1., 1., 3., 5., 7., 10.
239   };
240   const Int_t nVertexBins(sizeof(vertexBins)/sizeof(Double_t) - 1);
241
242   fPoolMgr = new AliEventPoolManager(poolsize, trackDepth,
243                                      nCentralityBins, centralityBins,
244                                      nVertexBins,     vertexBins);
245 }
246
247 THnSparse* AliAnalysisTaskLongRangeCorrelations::MakeHistSparsePhiEta(const char* name) const {
248   const Int_t   nBinsM1[] = { fnBinsPhi, fnBinsEta };
249   const Double_t xMinM1[] = {  fxMinPhi,  fxMinEta };
250   const Double_t xMaxM1[] = {  fxMaxPhi,  fxMaxEta };
251   const TString title(TString(name)
252                       +";#phi;#eta;");
253   return new THnSparseD(name, title.Data(), 2, nBinsM1, xMinM1, xMaxM1);
254 }
255 THnSparse* AliAnalysisTaskLongRangeCorrelations::MakeHistSparsePhiEtaPhiEta(const char* name) const {
256   const Int_t   nBinsM2[] = {  fnBinsPhi, fnBinsEta,  fnBinsPhi, fnBinsEta };
257   const Double_t xMinM2[] = {  fxMinPhi,  fxMinEta,   fxMinPhi,  fxMinEta  };
258   const Double_t xMaxM2[] = {  fxMaxPhi,  fxMaxEta,   fxMaxPhi,  fxMaxEta  };
259   const TString title(TString(name)
260                       +";#phi_{1};#eta_{1}"
261                       +";#phi_{2};#eta_{2};");
262   return new THnSparseD(name, title.Data(), 4, nBinsM2, xMinM2, xMaxM2);
263 }
264
265 TObjArray* AliAnalysisTaskLongRangeCorrelations::GetAcceptedTracks(AliAODEvent* pAOD,
266                                                                    Double_t centrality) {
267   TObjArray* tracks= new TObjArray;
268   tracks->SetOwner(kTRUE);
269
270   const Long64_t N(pAOD->GetNumberOfTracks());
271   AliDebug(5, Form("#tracks= %6lld %f", N, centrality));
272   for (Long64_t i(0); i<N; ++i) {
273     AliAODTrack* pAODTrack(dynamic_cast<AliAODTrack*>(pAOD->GetTrack(i)));
274     if (NULL == pAODTrack) continue;
275
276     // track filter selection
277     if (!pAODTrack->TestFilterBit(fTrackFilter)) continue;
278
279     // select only primary tracks
280     if (pAODTrack->GetType() != AliAODTrack::kPrimary) continue;
281
282     // select only charged tracks
283     if (pAODTrack->Charge() == 0) continue;
284
285     Fill("histQACentPt", pAODTrack->Charge()>0, centrality,       pAODTrack->Pt());
286     Fill("histQAPhiEta", pAODTrack->Charge()>0, pAODTrack->Phi(), pAODTrack->Eta());
287     if (pAODTrack->Phi() < fPhiMin || pAODTrack->Phi() > fPhiMax) continue;
288     if (pAODTrack->Pt()  < fPtMin  || pAODTrack->Pt()  > fPtMax)  continue;
289
290     tracks->Add(new LRCParticle(pAODTrack->Eta(), pAODTrack->Phi()));
291   } // next track
292   return tracks;
293 }
294
295 TObjArray* AliAnalysisTaskLongRangeCorrelations::GetAcceptedTracks(TClonesArray* tracksMC,
296                                                                    Double_t centrality) {
297   TObjArray* tracks= new TObjArray;
298   tracks->SetOwner(kTRUE);
299
300   const Long64_t N(tracksMC->GetEntriesFast());
301   AliDebug(5, Form("#tracks= %6lld %f", N, centrality));
302   for (Long64_t i(0); i<N; ++i) {
303     AliAODMCParticle* pMCTrack(dynamic_cast<AliAODMCParticle*>(tracksMC->At(i)));
304     if (NULL == pMCTrack) continue;    
305
306     // no track filter selection for MC tracks    
307
308     // select only primary tracks
309     if (kFALSE == pMCTrack->IsPhysicalPrimary()) continue;
310
311     // select only charged tracks
312     if (pMCTrack->Charge() == 0) continue;
313
314     Fill("MC_histQACentPt", pMCTrack->Charge()>0, centrality,      pMCTrack->Pt());
315     Fill("MC_histQAPhiEta", pMCTrack->Charge()>0, pMCTrack->Phi(), pMCTrack->Eta());
316
317     if (pMCTrack->Phi() < fPhiMin || pMCTrack->Phi() > fPhiMax) continue;
318     if (pMCTrack->Pt()  < fPtMin  || pMCTrack->Pt()  > fPtMax)  continue;
319
320     tracks->Add(new LRCParticle(pMCTrack->Eta(), pMCTrack->Phi()));
321   } // next track
322   return tracks;
323 }
324
325 void AliAnalysisTaskLongRangeCorrelations::CalculateMoments(TString prefix,
326                                                             TObjArray* tracks1,
327                                                             TObjArray* tracks2,
328                                                             Double_t weight) {
329   THnSparse* hN1(dynamic_cast<THnSparse*>(fOutputList->FindObject(prefix+"histMoment1PhiEta_1")));
330   THnSparse* hN2(dynamic_cast<THnSparse*>(fOutputList->FindObject(prefix+"histMoment1PhiEta_2")));
331   if (NULL == hN1) return;
332   if (NULL == hN2) return;
333
334   // <n_1>
335   THnSparse* hN1ForThisEvent(ComputeNForThisEvent(tracks1, "hN1"));
336   hN1->Add(hN1ForThisEvent, weight);
337
338   // <n_2>
339   THnSparse* hN2ForThisEvent(ComputeNForThisEvent(tracks2, "hN2"));
340   hN2->Add(hN2ForThisEvent, weight);
341
342   // n(eta) distributions
343   FillNEtaHist(prefix+"histNEta_300", hN1ForThisEvent, weight);
344   FillNEtaHist(prefix+"histNEta_120", hN1ForThisEvent, weight);
345   FillNEtaHist(prefix+"histNEta__30", hN1ForThisEvent, weight);
346   FillNEtaHist(prefix+"histNEta__15", hN1ForThisEvent, weight);
347
348   TObjArray* hNs(new TObjArray);
349
350   // <n_1 n_1>
351   hNs->AddAt(hN1ForThisEvent, 0);
352   hNs->AddAt(hN1ForThisEvent, 1);
353   ComputeNXForThisEvent(hNs, prefix+"histMoment2PhiEtaPhiEta_11", weight);
354
355   // <n_1 n_2>
356   hNs->AddAt(hN2ForThisEvent, 1);
357   ComputeNXForThisEvent(hNs, prefix+"histMoment2PhiEtaPhiEta_12", weight);
358
359   // <n_2 n_2>
360   hNs->AddAt(hN2ForThisEvent, 0);
361   ComputeNXForThisEvent(hNs, prefix+"histMoment2PhiEtaPhiEta_22", weight);
362
363   // clean up
364   delete hNs;
365   delete hN1ForThisEvent;
366   delete hN2ForThisEvent;
367 }
368
369 void AliAnalysisTaskLongRangeCorrelations::FillNEtaHist(TString name,
370                                                         THnSparse* hs,
371                                                         Double_t weight) {
372
373   TH2* hSum(dynamic_cast<TH2*>(fOutputList->FindObject(name)));
374   if (NULL == hSum) return;
375
376   TH2* hPerEvent(dynamic_cast<TH2*>(hSum->Clone("hPerEvent")));
377   if (NULL == hPerEvent) return;
378   hPerEvent->Reset();  
379
380   // fill hPerEvent
381   const Long64_t N(hs->GetNbins());
382   for (Long64_t i(0); i<N; ++i) {
383     Int_t coord[2] = { 0, 0 };
384     const Double_t n(hs->GetBinContent(i, coord));
385     const Double_t eta(hs->GetAxis(1)->GetBinCenter(coord[1]));
386     hPerEvent->Fill(eta, n);
387   }
388
389   // add zero counts for those eta bins with zero tracks
390   TH1* h(hPerEvent->ProjectionX());
391   for (Int_t i(1); i<=h->GetNbinsX(); ++i)
392     hPerEvent->SetBinContent(i,1, Double_t(h->GetBinContent(i) == 0));
393
394   hSum->Add(hPerEvent, weight);
395
396   delete h;
397   delete hPerEvent;
398 }
399
400 THnSparse* AliAnalysisTaskLongRangeCorrelations::ComputeNForThisEvent(TObjArray* tracks, const char* histName) const {
401   THnSparse* hN(MakeHistSparsePhiEta(histName));
402   const Long64_t nTracks(tracks->GetEntriesFast());
403   for (Long64_t i(0); i<nTracks; ++i) {
404     const LRCParticle* p(dynamic_cast<LRCParticle*>(tracks->At(i)));
405     if (NULL == p) continue;
406     const Double_t x[] = { p->Phi(), p->Eta() };
407     hN->Fill(x);
408   }
409   return hN;
410 }
411
412 class MultiDimIterator {
413 public:
414   MultiDimIterator(TObjArray* _fHs)
415     : fHs(_fHs)
416     , fN(fHs->GetEntriesFast())
417     , fDims(fN, 0)
418     , fIdxs(fN, 0)
419     , fNs  (fN, 0)
420     , fX (2*fN, 0)
421     , fj(0) {
422     for (Long64_t i=0; i<fN; ++i) {
423       THnSparse* hNi(reinterpret_cast<THnSparse*>(fHs->At(i)));
424       fDims[i] = hNi->GetNbins();
425       update(i);
426     }
427   }
428
429   const Double_t* GetX() const { return &fX.front(); }
430   Double_t        GetN() const {
431     return std::accumulate(fNs.begin(), fNs.end(), Double_t(1), std::multiplies<Double_t>());
432   }
433   Bool_t end() const { return fj == fN; }
434
435   MultiDimIterator& operator++() {
436     Long64_t j(0);
437     for (; j<fN; ++j) {
438       ++fIdxs[j];
439       update(j);
440       if (fIdxs[j] < fDims[j]) break;
441       fIdxs[j] = 0;
442       update(j);
443     }
444     fj = j;
445     return *this;
446   }
447
448 protected:
449   void update(Long64_t j) {
450     fj = j;
451     THnSparse* hs(reinterpret_cast<THnSparse*>(fHs->At(j)));
452     Int_t coord[] = { 0, 0 };
453     fNs[j] = hs->GetBinContent(fIdxs[j], coord);
454     for (Int_t k(0); k<2; ++k)
455       fX[2*j+k] = hs->GetAxis(k)->GetBinCenter(coord[k]);
456   }
457 private:
458   MultiDimIterator(const MultiDimIterator&);
459   MultiDimIterator& operator=(const MultiDimIterator&);
460
461   TObjArray*            fHs;
462   const Long64_t        fN;
463   std::vector<Long64_t> fDims;
464   std::vector<Long64_t> fIdxs;
465   std::vector<Double_t> fNs;
466   std::vector<Double_t> fX;
467   Long64_t              fj;
468 } ;
469
470 void AliAnalysisTaskLongRangeCorrelations::ComputeNXForThisEvent(TObjArray* hNs,
471                                                                  const char* histName,
472                                                                  Double_t weight) {
473   if (NULL == fOutputList) return;
474
475   THnSparse* hs(dynamic_cast<THnSparse*>(fOutputList->FindObject(histName)));
476   if (hs == NULL) return;
477
478   for (MultiDimIterator mdi(hNs); !mdi.end(); ++mdi)
479     hs->Fill(mdi.GetX(), mdi.GetN()*weight);
480 }
481
482 void AliAnalysisTaskLongRangeCorrelations::Fill(const char* histName, Double_t x) {
483   if (NULL == fOutputList) return;
484   TH1* h = dynamic_cast<TH1*>(fOutputList->FindObject(histName));
485   if (h == NULL) return;
486   h->Fill(x);
487 }
488 void AliAnalysisTaskLongRangeCorrelations::Fill(const char* histName, Double_t x, Double_t y) {
489   if (NULL == fOutputList) return;
490   TH2* h = dynamic_cast<TH2*>(fOutputList->FindObject(histName));
491   if (h == NULL) return;
492   h->Fill(x, y);
493 }
494 void AliAnalysisTaskLongRangeCorrelations::Fill(const char* histName, Double_t x, Double_t y, Double_t z) {
495   if (NULL == fOutputList) return;
496   TH3* h = dynamic_cast<TH3*>(fOutputList->FindObject(histName));
497   if (h == NULL) return;
498   h->Fill(x, y, z);
499 }
500 void AliAnalysisTaskLongRangeCorrelations::Fill(const char* histName, const Double_t* x, Double_t w) {
501   if (NULL == fOutputList) return;
502   THnSparse* h = dynamic_cast<THnSparse*>(fOutputList->FindObject(histName));
503   if (h == NULL) return;
504   h->Fill(x, w);
505 }