]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/CaloCalib/AliAnalysisTaskEMCALPi0PbPb.cxx
bugfix for event counter
[u/mrichter/AliRoot.git] / PWG4 / CaloCalib / AliAnalysisTaskEMCALPi0PbPb.cxx
1 // $Id$
2
3 #include "AliAnalysisTaskEMCALPi0PbPb.h"
4 #include <TAxis.h>
5 #include <TChain.h>
6 #include <TClonesArray.h>
7 #include <TFile.h>
8 #include <TH1F.h>
9 #include <TH2F.h>
10 #include <TList.h>
11 #include <TLorentzVector.h>
12 #include <TNtuple.h>
13 #include "AliAODEvent.h"
14 #include "AliAODVertex.h"
15 #include "AliAnalysisManager.h"
16 #include "AliAnalysisTaskEMCALClusterizeFast.h"
17 #include "AliCentrality.h"
18 #include "AliEMCALGeoUtils.h"
19 #include "AliESDEvent.h"
20 #include "AliESDVertex.h"
21 #include "AliLog.h"
22
23 ClassImp(AliAnalysisTaskEMCALPi0PbPb)
24
25 //________________________________________________________________________
26 AliAnalysisTaskEMCALPi0PbPb::AliAnalysisTaskEMCALPi0PbPb(const char *name) 
27   : AliAnalysisTaskSE(name),
28     fAsymMax(1),
29     fCentVar("V0M"),
30     fCentFrom(0),
31     fCentTo(100),
32     fVtxZMin(-10),
33     fVtxZMax(+10),
34     fUseQualFlag(1),
35     fClusName(),
36     fDoNtuple(0),
37     fDoAfterburner(0),
38     fNminCells(1),
39     fGeom(0),
40     fOutput(0),
41     fEsdEv(0),
42     fAodEv(0),
43     fRecPoints(0),
44     fEsdClusters(0),
45     fEsdCells(0),
46     fAodClusters(0),
47     fAodCells(0),
48     fPtRanges(0),
49     fNtuple(0),
50     fHCuts(0x0),
51     fHVertexZ(0x0),
52     fHVertexZ2(0x0),
53     fHCent(0x0),
54     fHCentQual(0x0),
55     fHColuRow(0x0),
56     fHColuRowE(0x0),
57     fHCellMult(0x0),
58     fHCellE(0x0),
59     fHCellH(0x0),
60     fHCellM(0x0),
61     fHCellM2(0x0),
62     fHClustEccentricity(0),
63     fHClustEtaPhi(0x0),
64     fHClustEnergyPt(0x0),
65     fHClustEnergySigma(0x0),
66     fHClustSigmaSigma(0x0),
67     fHClustNCellEnergyRatio(0x0),
68     fHPionEtaPhi(0x0),
69     fHPionMggPt(0x0),
70     fHPionMggAsym(0x0),
71     fHPionMggDgg(0x0)
72 {
73   // Constructor.
74
75   if (!name)
76     return;
77   SetName(name);
78   DefineInput(0, TChain::Class());
79   DefineOutput(1, TList::Class());
80   fBranchNames="ESD:AliESDRun.,AliESDHeader.,PrimaryVertex.,SPDVertex.,TPCVertex.,EMCALCells. "
81                "AOD:header,vertices,emcalCells";
82 }
83
84 //________________________________________________________________________
85 AliAnalysisTaskEMCALPi0PbPb::~AliAnalysisTaskEMCALPi0PbPb()
86 {
87   // Destructor.
88
89   delete fOutput; fOutput = 0;
90   delete fPtRanges; fPtRanges = 0;
91   delete fGeom; fGeom = 0;
92   delete [] fHColuRow;
93   delete [] fHColuRowE;
94   delete [] fHCellMult;
95 }
96
97 //________________________________________________________________________
98 void AliAnalysisTaskEMCALPi0PbPb::UserCreateOutputObjects()
99 {
100   // Create user objects here.
101
102   fOutput = new TList();
103   fOutput->SetOwner();
104
105   if (fDoNtuple) {
106     TFile *f = OpenFile(1);
107     if (f)
108       fNtuple = new TNtuple(Form("nt%.0fto%.0f",fCentFrom,fCentTo),"nt",
109                             "run:evt:cent:pt:e:emax:n:db:disp:mn:ms:chi:cpv:ecc:sig:eta:phi");
110   }
111
112   // histograms
113   TH1::SetDefaultSumw2(kTRUE);
114   TH2::SetDefaultSumw2(kTRUE);
115   fHCuts = new TH1F("hEventCuts","",4,0.5,4.5);
116   fHCuts->GetXaxis()->SetBinLabel(1,"All (PS)");
117   fHCuts->GetXaxis()->SetBinLabel(2,Form("%s: %.0f-%.0f",fCentVar.Data(),fCentFrom,fCentTo));
118   fHCuts->GetXaxis()->SetBinLabel(3,"QFlag");
119   fHCuts->GetXaxis()->SetBinLabel(4,Form("zvtx: %.0f-%.0f",fVtxZMin,fVtxZMax));
120   fOutput->Add(fHCuts);
121   fHVertexZ = new TH1F("hVertexZBeforeCut","",100,-25,25);
122   fHVertexZ->SetXTitle("z [cm]");
123   fOutput->Add(fHVertexZ);
124   fHVertexZ2 = new TH1F("hVertexZAfterCut","",100,-25,25);
125   fHVertexZ2->SetXTitle("z [cm]");
126   fOutput->Add(fHVertexZ2);
127   fHCent = new TH1F("hCentBeforeCut","",101,-1,100);
128   fHCent->SetXTitle(fCentVar.Data());
129   fOutput->Add(fHCent);
130   fHCentQual = new TH1F("hCentAfterCut","",101,-1,100);
131   fHCentQual->SetXTitle(fCentVar.Data());
132   fOutput->Add(fHCentQual);
133
134   // histograms for cells
135   fHColuRow   = new TH2F*[4];
136   fHColuRowE  = new TH2F*[4];
137   fHCellMult  = new TH1F*[4];
138   for (Int_t i = 0; i<4; ++i) {
139     fHColuRow[i] = new TH2F(Form("hColRow_Mod%d",i),"",49,0,49,25,0,25);
140     fHColuRow[i]->SetTitle(Form("Module %d: Occupancy", i));
141     fHColuRow[i]->SetXTitle("col (i#eta)");
142     fHColuRow[i]->SetYTitle("row (i#phi)");
143     fHColuRowE[i] = new TH2F(Form("hColRowE_Mod%d", i),"",49,0,49,25,0,25);
144     fHColuRowE[i]->SetTitle(Form("Module %d: Cell energy",i));
145     fHColuRowE[i]->SetXTitle("col (i#eta)");
146     fHColuRowE[i]->SetYTitle("row (i#phi)");
147     fHCellMult[i] = new TH1F(Form("hCellMult_Mod%d",i),"",1000,0,1000); 
148     fHCellMult[i]->SetTitle(Form("Module %d: Cell multiplicity",i));
149     fHCellMult[i]->SetXTitle("# of cells");
150     fOutput->Add(fHColuRow[i]);
151     fOutput->Add(fHColuRowE[i]);
152     fOutput->Add(fHCellMult[i]);
153   }  
154   fHCellE = new TH1F("hCellE","",150,0.,15.);
155   fHCellE->SetXTitle("E_{cell} [GeV]");
156   fOutput->Add(fHCellE);
157   fHCellH = new TH1F ("fHCellHighestE","",150,0.,15.);
158   fHCellH->SetXTitle("E^{max}_{cell} [GeV]");
159   fOutput->Add(fHCellH);
160   fHCellM = new TH1F ("fHCellMeanEperHitCell","",250,0.,2.5);
161   fHCellM->SetXTitle("#LT E_{cell}#GT [GeV]");
162   fOutput->Add(fHCellM);
163   fHCellM2 = new TH1F ("fHCellMeanEperAllCells","",250,0.,1);
164   fHCellM2->SetXTitle("1/N_{cells} #Sum E_{cell} [GeV]");
165   fOutput->Add(fHCellM2);
166
167   // histograms for clusters
168   fHClustEccentricity = new TH1F("hClustEccentricity","",100,-0.1,1.1);
169   fHClustEccentricity->SetXTitle("#epsilon_{C}");
170   fOutput->Add(fHClustEccentricity);
171   fHClustEtaPhi = new TH2F("hClustEtaPhi","",500,-0.8,0.8,500,1.2,2.2);
172   fHClustEtaPhi->SetXTitle("#eta");
173   fHClustEtaPhi->SetYTitle("#varphi");
174   fOutput->Add(fHClustEtaPhi);
175   fHClustEnergyPt = new TH2F("hClustEnergyPt","",250,0,50,250,0,50);
176   fHClustEnergyPt->SetXTitle("E [GeV]");
177   fHClustEnergyPt->SetYTitle("p_{T} [GeV/c]");
178   fOutput->Add(fHClustEnergyPt);
179   fHClustEnergySigma = new TH2F("hClustEnergySigma","",100,0,100,500,0,50);
180   fHClustEnergySigma->SetXTitle("E_{C} * #sigma_{max} [GeV*cm]");
181   fHClustEnergySigma->SetYTitle("E_{C} [GeV]");
182   fOutput->Add(fHClustEnergySigma);
183   fHClustSigmaSigma = new TH2F("hClustSigmaSigma","",500,0,50,500,0,50);
184   fHClustSigmaSigma->SetXTitle("#lambda_{0} [cm]");
185   fHClustSigmaSigma->SetYTitle("#sigma_{max} [cm]");
186   fOutput->Add(fHClustSigmaSigma);
187   fHClustNCellEnergyRatio = new TH2F("hClustNCellEnergyRatio","",27,-0.5,26.5,101,-0.05,1.05);
188   fHClustNCellEnergyRatio->SetXTitle("N_{cells}");
189   fHClustNCellEnergyRatio->SetYTitle("E^{max}_{cell}/E_{clus}");
190   fOutput->Add(fHClustNCellEnergyRatio);
191
192   // histograms for pion candidates
193   fHPionEtaPhi = new TH2F("hPionEtaPhi","",100,-0.8,0.8,100,1.2,2.2);
194   fHPionEtaPhi->SetXTitle("#eta_{#gamma#gamma}");
195   fHPionEtaPhi->SetYTitle("#varphi_{#gamma#gamma}");
196   fOutput->Add(fHPionEtaPhi);
197   fHPionMggPt = new TH2F("hPionMggPt","",1000,0,2,100,0,20.0);
198   fHPionMggPt->SetXTitle("M_{#gamma#gamma} [GeV/c^{2}]");
199   fHPionMggPt->SetYTitle("p_{T}^{#gamma#gamma} [GeV/c]");
200   fOutput->Add(fHPionMggPt);
201   fHPionMggAsym = new TH2F("hPionMggAsym","",1000,0,2,100,0,1);
202   fHPionMggAsym->SetXTitle("M_{#gamma#gamma} [GeV/c^{2}]");
203   fHPionMggAsym->SetYTitle("Z_{#gamma#gamma} [GeV]");
204   fOutput->Add(fHPionMggAsym);
205   fHPionMggDgg = new TH2F("hPionMggDgg","",1000,0,2,100,0,10);
206   fHPionMggDgg->SetXTitle("M_{#gamma#gamma} [GeV/c^{2}]");
207   fHPionMggDgg->SetYTitle("opening angle [grad]");
208   fOutput->Add(fHPionMggDgg);
209   const Int_t nbins = 19; 
210   Double_t xbins[nbins] = {0.5,1,1.5,2,2.5,3,3.5,4,4.5,6,7,8,9,10,12.5,15,20,25,50};
211   fPtRanges = new TAxis(nbins-1,xbins);
212   for (Int_t i = 0; i<=nbins; ++i) {
213     fHPionInvMasses[i] = new TH1F(Form("hPionInvMass%d",i),"",1000,0,2);
214     fHPionInvMasses[i]->SetXTitle("M_{#gamma#gamma} [GeV/c^{2}]");
215     if (i==0)
216       fHPionInvMasses[i]->SetTitle(Form("0 < p_{T}^{#gamma#gamma} <%.1f",xbins[0]));
217     else if (i==nbins)
218       fHPionInvMasses[i]->SetTitle(Form("p_{T}^{#gamma#gamma} > 50"));
219     else 
220       fHPionInvMasses[i]->SetTitle(Form("%.1f < p_{T}^{#gamma#gamma} <%.1f",xbins[i-1],xbins[i]));
221     fOutput->Add(fHPionInvMasses[i]);
222   }
223  
224   PostData(1, fOutput); 
225 }
226
227 //________________________________________________________________________
228 void AliAnalysisTaskEMCALPi0PbPb::UserExec(Option_t *) 
229 {
230   // Called for each event.
231
232   if (!InputEvent())
233     return;
234
235   AliAnalysisManager *am = AliAnalysisManager::GetAnalysisManager();
236   fEsdEv = dynamic_cast<AliESDEvent*>(InputEvent());
237   if (fEsdEv) {
238     am->LoadBranch("AliESDRun.");
239     am->LoadBranch("AliESDHeader.");
240   } else {
241     fAodEv = dynamic_cast<AliAODEvent*>(InputEvent());
242     am->LoadBranch("header");
243   }
244
245   if (!fGeom) { // set misalignment matrices (stored in first event)
246     fGeom = new AliEMCALGeoUtils("EMCAL_FIRSTYEARV1","EMCAL");
247     if (fEsdEv) {
248       for (Int_t i=0; i<4; ++i)
249         fGeom->SetMisalMatrix(fEsdEv->GetESDRun()->GetEMCALMatrix(i),i);
250     } else {
251       for (Int_t i=0; i<4; ++i)
252         fGeom->SetMisalMatrix(fAodEv->GetHeader()->GetEMCALMatrix(i),i);
253     }
254     fGeom->GetEMCGeometry();
255   }
256
257   Int_t cut = 1;
258   fHCuts->Fill(cut++);
259
260   const AliCentrality *centP = InputEvent()->GetCentrality();
261   Double_t cent = centP->GetCentralityPercentileUnchecked(fCentVar);
262   fHCent->Fill(cent);
263   if (cent<fCentFrom||cent>fCentTo)
264     return;
265
266   fHCuts->Fill(cut++);
267   
268   if (fUseQualFlag) {
269     if (centP->GetQuality()>0)
270       return;
271   }
272
273   fHCentQual->Fill(cent);
274   fHCuts->Fill(cut++);
275
276   if (fEsdEv) {
277     am->LoadBranch("PrimaryVertex.");
278     am->LoadBranch("SPDVertex.");
279     am->LoadBranch("TPCVertex.");
280   } else {
281     fAodEv = dynamic_cast<AliAODEvent*>(InputEvent());
282     am->LoadBranch("vertices");
283   }
284
285   const AliVVertex *vertex = InputEvent()->GetPrimaryVertex();
286   if (!vertex)
287     return;
288
289   fHVertexZ->Fill(vertex->GetZ());
290
291   if(vertex->GetZ()<fVtxZMin||vertex->GetZ()>fVtxZMax)
292     return;
293
294   fHCuts->Fill(cut++);
295   fHVertexZ2->Fill(vertex->GetZ());
296
297   fRecPoints   = 0; // will be set if fClusName is given and AliAnalysisTaskEMCALClusterizeFast is used
298   fEsdClusters = 0; // will be set if ESD input used and if fRecPoints are not set of if clusters are attached
299   fEsdCells    = 0; // will be set if ESD input used
300   fAodClusters = 0; // will be set if AOD input used and if fRecPoints are not set of if clusters are attached
301                     //             or if fClusName is given and AliAnalysisTaskEMCALClusterizeFast in AOD output mode
302   fAodCells    = 0; // will be set if AOD input used
303
304   // deal with special output from AliAnalysisTaskEMCALClusterizeFast first
305   Bool_t clusattached = 0;
306   Bool_t recalibrated = 0;
307   if (1 && !fClusName.IsNull()) {
308     AliAnalysisTaskEMCALClusterizeFast *cltask = 0;
309     TObjArray *ts = am->GetTasks();
310     cltask = dynamic_cast<AliAnalysisTaskEMCALClusterizeFast*>(ts->FindObject(fClusName));
311     if (cltask && cltask->GetClusters()) {
312       fRecPoints = const_cast<TObjArray*>(cltask->GetClusters());
313       clusattached = cltask->GetAttachClusters();
314       if (cltask->GetCalibData()!=0)
315         recalibrated = kTRUE;
316     }
317   }
318   if (1 && AODEvent() && !fClusName.IsNull()) {
319     TList *l = AODEvent()->GetList();
320     TClonesArray *clus = 0;
321     if (l) {
322       clus = dynamic_cast<TClonesArray*>(l->FindObject(fClusName));
323       fAodClusters = clus;
324     }
325   }
326
327   if (fEsdEv) { // ESD input mode
328     if (1 && (!fRecPoints||clusattached)) {
329       if (!clusattached)
330         am->LoadBranch("CaloClusters");
331       TList *l = fEsdEv->GetList();
332       TClonesArray *clus = 0;
333       if (l) {
334         clus = dynamic_cast<TClonesArray*>(l->FindObject("CaloClusters"));
335         fEsdClusters = clus;
336       }
337     }
338     if (1) {
339       if (!recalibrated)
340         am->LoadBranch("EMCALCells.");
341       fEsdCells = fEsdEv->GetEMCALCells();
342     }
343   } else if (fAodEv) { // AOD input mode
344     if (1 && (!fAodClusters || clusattached)) {
345       if (!clusattached)
346         am->LoadBranch("caloClusters");
347       TList *l = fAodEv->GetList();
348       TClonesArray *clus = 0;
349       if (l) {
350         clus = dynamic_cast<TClonesArray*>(l->FindObject("caloClusters"));
351         fAodClusters = clus;
352       }
353     }
354     if (1) {
355       if (!recalibrated)
356         am->LoadBranch("emcalCells");
357       fAodCells = fAodEv->GetEMCALCells();
358     }
359   } else {
360     AliFatal("Impossible to not have either pointer to ESD or AOD event");
361   }
362
363   if (1) {
364     AliDebug(2,Form("fRecPoints   set: %p", fRecPoints));
365     AliDebug(2,Form("fEsdClusters set: %p", fEsdClusters));
366     AliDebug(2,Form("fEsdCells    set: %p", fEsdCells));
367     AliDebug(2,Form("fAodClusters set: %p", fAodClusters));
368     AliDebug(2,Form("fAodCells    set: %p", fAodCells));
369   }
370
371   if (fDoAfterburner)
372     ClusterAfterburner();
373
374   FillCellHists();
375   FillClusHists();
376   FillPionHists();
377
378   PostData(1, fOutput);
379 }      
380
381 //________________________________________________________________________
382 void AliAnalysisTaskEMCALPi0PbPb::Terminate(Option_t *) 
383 {
384   // Terminate called at the end of analysis.
385
386   if (fNtuple) {
387     TFile *f = OpenFile(1);
388     if (f) 
389       fNtuple->Write();
390   }
391 }
392
393 //________________________________________________________________________
394 void AliAnalysisTaskEMCALPi0PbPb::FillCellHists()
395 {
396   // Fill histograms related to cell properties.
397   
398   AliVCaloCells *cells = fEsdCells;
399   if (!cells)
400     cells = fAodCells;
401
402   if (cells) {
403     Int_t cellModCount[4] = {0,0,0,0};
404     Double_t cellMaxE = 0; 
405     Double_t cellMeanE = 0; 
406     Int_t ncells = cells->GetNumberOfCells();
407     if (ncells==0)
408       return;
409
410     for (Int_t i = 0; i<ncells; ++i ) {
411       Int_t absID    = cells->GetCellNumber(i);
412       Double_t cellE = cells->GetAmplitude(i);
413       fHCellE->Fill(cellE);
414       if (cellE>cellMaxE) 
415         cellMaxE = cellE;
416       cellMeanE += cellE;
417
418       Int_t iSM=-1, iTower=-1, nIphi=-1, nIeta=-1;
419       Bool_t ret = fGeom->GetCellIndex(absID, iSM, iTower, nIphi, nIeta);
420       if (!ret) {
421         AliError(Form("Could not get cell index for %d", absID));
422         continue;
423       }
424       ++cellModCount[iSM];
425       Int_t iPhi=-1, iEta=-1;
426       fGeom->GetCellPhiEtaIndexInSModule(iSM, iTower, nIphi, nIeta, iPhi, iEta);   
427       fHColuRow[iSM]->Fill(iEta,iPhi,1);
428       fHColuRowE[iSM]->Fill(iEta,iPhi,cellE);
429     }    
430     fHCellH->Fill(cellMaxE);
431     cellMeanE /= ncells;
432     fHCellM->Fill(cellMeanE);
433     fHCellM2->Fill(cellMeanE*ncells/24/48/4); //hard-coded but there is a way to figure out from geometry
434     for (Int_t i=0; i<4; ++i) 
435       fHCellMult[i]->Fill(cellModCount[i]);
436   }
437 }
438
439 //________________________________________________________________________
440 void AliAnalysisTaskEMCALPi0PbPb::FillClusHists()
441 {
442   // Fill histograms related to cluster properties.
443   Double_t clusterEcc = 0;
444
445
446   Double_t vertex[3] = {0,0,0};
447   InputEvent()->GetPrimaryVertex()->GetXYZ(vertex);
448
449   TObjArray *clusters = fEsdClusters;
450   if (!clusters)
451     clusters = fAodClusters;
452
453   if (clusters) {
454     TLorentzVector clusterVec;
455     Int_t nclus = clusters->GetEntries();
456     for(Int_t i = 0; i<nclus; ++i) {
457       AliVCluster *clus = static_cast<AliVCluster*>(clusters->At(i));
458       if (!clus)
459         continue;
460       if (!clus->IsEMCAL()) 
461         continue;
462       clus->GetMomentum(clusterVec,vertex);
463       Double_t maxAxis = 0, minAxis = 0;
464       GetSigma(clus,maxAxis,minAxis);
465       if (maxAxis > 0)
466         clusterEcc = TMath::Sqrt(1.0 - minAxis*minAxis/(maxAxis*maxAxis));
467
468       fHClustEccentricity->Fill(clusterEcc); 
469       fHClustEtaPhi->Fill(clusterVec.Eta(),clusterVec.Phi());
470       fHClustEnergyPt->Fill(clusterVec.E(),clusterVec.Pt());
471       fHClustEnergySigma->Fill(clus->E()*maxAxis,clus->E());
472       fHClustSigmaSigma->Fill(max(clus->GetM02(),clus->GetM20()),clus->E()*maxAxis);
473       fHClustNCellEnergyRatio->Fill(clus->GetNCells(),GetMaxCellEnergy(clus)/clus->E());
474       if (fNtuple) {
475         Float_t vals[18];
476         vals[0]  = InputEvent()->GetRunNumber();
477         vals[1]  = (((UInt_t)InputEvent()->GetOrbitNumber()  << 12) | (UInt_t)InputEvent()->GetBunchCrossNumber()); 
478         vals[2]  = InputEvent()->GetCentrality()->GetCentralityPercentileUnchecked(fCentVar);
479         vals[3]  = clusterVec.Pt();
480         vals[4]  = clusterVec.E();
481         vals[5]  = GetMaxCellEnergy(clus);
482         vals[6]  = clus->GetNCells();
483         vals[7]  = clus->GetDistanceToBadChannel();
484         vals[8]  = clus->GetDispersion();
485         vals[9]  = clus->GetM20();
486         vals[10] = clus->GetM02();
487         vals[11] = clus->Chi2();
488         vals[12] = clus->GetEmcCpvDistance();
489         vals[13] = clusterEcc;
490         vals[14] = GetMaxCellEnergy(clus)/clus->E();
491         vals[15] = maxAxis;
492         vals[16] = clusterVec.Eta();
493         vals[17] = clusterVec.Phi();
494         fNtuple->Fill(vals);
495       }
496     }
497   }
498 }
499
500 //________________________________________________________________________
501 void AliAnalysisTaskEMCALPi0PbPb::FillPionHists()
502 {
503   // Fill histograms related to pions.
504
505   Double_t vertex[3] = {0,0,0};
506   InputEvent()->GetPrimaryVertex()->GetXYZ(vertex);
507
508   TObjArray *clusters = fEsdClusters;
509   if (!clusters)
510     clusters = fAodClusters;
511
512   if (clusters) {
513     TLorentzVector clusterVec1;
514     TLorentzVector clusterVec2;
515     TLorentzVector pionVec;
516
517     Int_t nclus = clusters->GetEntries();
518     for (Int_t i = 0; i<nclus; ++i) {
519       AliVCluster *clus1 = static_cast<AliVCluster*>(clusters->At(i));
520       if (!clus1)
521         continue;
522       if (!clus1->IsEMCAL()) 
523         continue;
524       if (clus1->E()<0.010)
525         continue;
526       if (clus1->GetNCells()<fNminCells)
527         continue;
528       clus1->GetMomentum(clusterVec1,vertex);
529       for (Int_t j = i+1; j<nclus; ++j) {
530         AliVCluster *clus2 = static_cast<AliVCluster*>(clusters->At(j));
531         if (!clus2)
532           continue;
533         if (!clus2->IsEMCAL()) 
534           continue;
535         if (clus2->E()<0.010)
536           continue;
537         if (clus2->GetNCells()<fNminCells)
538           continue;
539         clus2->GetMomentum(clusterVec2,vertex);
540         pionVec = clusterVec1 + clusterVec2;
541         Double_t pionZgg = TMath::Abs(clusterVec1.E()-clusterVec2.E())/pionVec.E();
542         Double_t pionOpeningAngle = clusterVec1.Angle(clusterVec2.Vect());
543         if (pionZgg < fAsymMax) {
544           fHPionEtaPhi->Fill(pionVec.Eta(),pionVec.Phi()); 
545           fHPionMggPt->Fill(pionVec.M(),pionVec.Pt()); 
546           fHPionMggAsym->Fill(pionVec.M(),pionZgg); 
547           fHPionMggDgg->Fill(pionVec.M(),pionOpeningAngle); 
548           Int_t bin = fPtRanges->FindBin(pionVec.Pt());
549           fHPionInvMasses[bin]->Fill(pionVec.M());
550         }
551       }
552     }
553   } 
554 }
555
556 //________________________________________________________________________
557 void  AliAnalysisTaskEMCALPi0PbPb::ClusterAfterburner()
558 {
559   // 
560
561   AliVCaloCells *cells = fEsdCells;
562   if (!cells)
563     cells = fAodCells;
564
565   if (!cells)
566     return;
567
568   Int_t ncells = cells->GetNumberOfCells();
569   if (ncells<=0)
570     return;
571
572   Double_t cellMeanE = 0, cellSigE = 0;
573   for (Int_t i = 0; i<ncells; ++i ) {
574     Double_t cellE = cells->GetAmplitude(i);
575     cellMeanE += cellE;
576     cellSigE += cellE*cellE;
577   }    
578   cellMeanE /= ncells;
579   cellSigE /= ncells;
580   cellSigE -= cellMeanE*cellMeanE;
581   if (cellSigE<0)
582     cellSigE = 0;
583   cellSigE = TMath::Sqrt(cellSigE / ncells);
584
585   Double_t subE = cellMeanE - 7*cellSigE;
586   if (subE<0)
587     return;
588
589   for (Int_t i = 0; i<ncells; ++i) {
590     Short_t id=-1;
591     Double_t amp=0,time=0;
592     if (!cells->GetCell(i, id, amp, time))
593       continue;
594     amp -= cellMeanE;
595     if (amp<0.001)
596       amp = 0;
597     cells->SetCell(i, id, amp, time);
598   }    
599
600   TObjArray *clusters = fEsdClusters;
601   if (!clusters)
602     clusters = fAodClusters;
603
604   if (!clusters)
605     return;
606
607   Int_t nclus = clusters->GetEntries();
608   for (Int_t i = 0; i<nclus; ++i) {
609     AliVCluster *clus = static_cast<AliVCluster*>(clusters->At(i));
610     Int_t nc = clus->GetNCells();
611     Double_t oldE = clus->E();
612     Double_t clusE = 0;
613     Short_t nc2 = 0;
614     UShort_t ids[100] = {0};
615     Double_t fra[100] = {0};
616     for (Int_t j = 0; j<nc; ++j) {
617       Short_t id = clus->GetCellAbsId(j);
618       Double_t cen = cells->GetCellAmplitude(id);
619       clusE += cen;
620       if (cen>0) {
621         ids[nc] = id;
622         ++nc;
623       }
624     }
625     if (clusE<=0) {
626       clusters->RemoveAt(i);
627       continue;
628     }
629     for (Int_t j = 0; j<nc2; ++j) {
630       Short_t id = ids[j];
631       Double_t cen = cells->GetCellAmplitude(id);
632       fra[j] = cen/clusE;
633     }
634     clus->SetE(clusE);
635     AliAODCaloCluster *aodclus = dynamic_cast<AliAODCaloCluster*>(clus);
636     if (aodclus) {
637       aodclus->Clear("");
638       aodclus->SetNCells(nc2);
639       aodclus->SetCellsAmplitudeFraction(fra);
640       aodclus->SetCellsAbsId(ids);
641     }
642     cout << i << " " << clusE << " " << oldE << endl;
643   }
644
645   clusters->Compress();
646 }
647
648 // virtual void     AddAt(TObject *obj, Int_t idx);
649 //   virtual TObject *RemoveAt(Int_t idx);
650
651 //________________________________________________________________________
652 Double_t  AliAnalysisTaskEMCALPi0PbPb::GetMaxCellEnergy(AliVCluster *cluster)
653 {
654   // Get maximum energy of attached cell.
655
656   Double_t maxe = 0;
657   Int_t ncells = cluster->GetNCells();
658   if (fEsdCells) {
659     for (Int_t i=0; i<ncells; i++) {
660       Double_t e = fEsdCells->GetCellAmplitude(cluster->GetCellAbsId(i));
661       if (e>maxe)
662         maxe = e;
663     }
664   } else {
665     for (Int_t i=0; i<ncells; i++) {
666       Double_t e = fAodCells->GetCellAmplitude(cluster->GetCellAbsId(i));
667       if (e>maxe)
668         maxe = e;
669     }
670   }
671   return maxe;
672 }
673
674 //________________________________________________________________________
675 void AliAnalysisTaskEMCALPi0PbPb::GetSigma(AliVCluster *cluster, Double_t& sigmaMax, Double_t &sigmaMin)
676 {
677   // Calculate the (E) weighted variance along the longer (eigen) axis.
678
679   sigmaMax = 0;               // cluster variance along its longer axis
680   sigmaMin = 0;               // cluster variance along its shorter axis
681   Double_t Ec = cluster->E(); // cluster energy
682   Double_t Xc = 0 ;           // cluster first moment along X
683   Double_t Yc = 0 ;           // cluster first moment along Y
684   Double_t Sxx = 0 ;          // cluster second central moment along X (variance_X^2)
685   Double_t Sxy = 0 ;          // cluster second central moment along Y (variance_Y^2)
686   Double_t Syy = 0 ;          // cluster covariance^2
687
688   AliVCaloCells *cells = fEsdCells;
689   if (!cells)
690     cells = fAodCells;
691
692   if (!cells)
693     return;
694
695   Int_t ncells = cluster->GetNCells();
696   if (ncells==1)
697     return;
698
699   TVector3 pos;
700   for(Int_t j=0; j<ncells; ++j) {
701     fGeom->GetGlobal(cluster->GetCellAbsId(j),pos);
702     Int_t id = cluster->GetCellAbsId(j);
703     Double_t cellen = cells->GetCellAmplitude(id);
704     Xc  += cellen*pos.X();
705     Yc  += cellen*pos.Y();    
706     Sxx += cellen*pos.X()*pos.X(); 
707     Syy += cellen*pos.Y()*pos.Y(); 
708     Sxy += cellen*pos.X()*pos.Y(); 
709   }
710   Xc  /= Ec;
711   Yc  /= Ec;
712   Sxx /= Ec;
713   Syy /= Ec;
714   Sxy /= Ec;
715   Sxx -= Xc*Xc;
716   Syy -= Yc*Yc;
717   Sxy -= Xc*Yc;
718   sigmaMax = (Sxx + Syy + TMath::Sqrt((Sxx-Syy)*(Sxx-Syy)+4.0*Sxy*Sxy))/2.0;
719   sigmaMax = TMath::Sqrt(sigmaMax); 
720   sigmaMin = TMath::Abs(Sxx + Syy - TMath::Sqrt((Sxx-Syy)*(Sxx-Syy)+4.0*Sxy*Sxy))/2.0;
721   sigmaMin = TMath::Sqrt(sigmaMin); 
722 }
723