]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGJE/EMCALJetTasks/AliAnalysisTaskScale.cxx
Remove hard-coded number of centrality bins
[u/mrichter/AliRoot.git] / PWGJE / EMCALJetTasks / AliAnalysisTaskScale.cxx
1 // $Id$
2 //
3 // Scale task.
4 //
5 // Author: R.Reed, M.Connors
6
7 #include "AliAnalysisTaskScale.h"
8
9 #include <TClonesArray.h>
10 #include <TF1.h>
11 #include <TH2F.h>
12 #include <TLorentzVector.h>
13 #include <TMath.h>
14
15 #include "AliEMCALGeometry.h"
16 #include "AliLog.h"
17 #include "AliVCluster.h"
18 #include "AliVParticle.h"
19 #include "AliParticleContainer.h"
20
21 ClassImp(AliAnalysisTaskScale)
22
23 //________________________________________________________________________
24 AliAnalysisTaskScale::AliAnalysisTaskScale() : 
25   AliAnalysisTaskEmcal("AliAnalysisTaskScale", kTRUE), 
26   fScaleFunction(0),
27   fEmcalArea(1),
28   fTpcArea(1),
29   fHistPtTPCvsCent(0), 
30   fHistPtEMCALvsCent(0), 
31   fHistEtvsCent(0),  
32   fHistScalevsCent(0),  
33   fHistDeltaScalevsCent(0), 
34   fHistScaleEmcalvsCent(0),      
35   fHistScale2EmcalvsCent(0),
36   fHistDeltaScale2EmcalvsCent(0),     
37   fHistChScalevsCent(0),          
38   fHistChScale2EmcalvsCent(0),   
39   fHistPtTPCvsNtrack(0), 
40   fHistPtEMCALvsNtrack(0), 
41   fHistEtvsNtrack(0),  
42   fHistScalevsNtrack(0),  
43   fHistDeltaScalevsNtrack(0),
44   fHistScaleEmcalvsNtrack(0),      
45   fHistScale2EmcalvsNtrack(0),     
46   fHistChScalevsNtrack(0),          
47   fHistChScale2EmcalvsNtrack(0),   
48   fHistTrackPtvsCent(0), 
49   fHistClusterPtvsCent(0),
50   fHistTrackEtaPhi(0), 
51   fHistClusterEtaPhi(0),
52   fHistScalevsScale2Emcal(0),      
53   fHistScalevsScaleEmcal(0),       
54   fHistScaleEmcalvsScale2Emcal(0)
55 {
56   // Default constructor.
57
58   SetMakeGeneralHistograms(kTRUE);
59 }
60
61 //________________________________________________________________________
62 AliAnalysisTaskScale::AliAnalysisTaskScale(const char *name) :
63   AliAnalysisTaskEmcal(name, kTRUE), 
64   fScaleFunction(0),
65   fEmcalArea(1),
66   fTpcArea(1),
67   fHistPtTPCvsCent(0), 
68   fHistPtEMCALvsCent(0), 
69   fHistEtvsCent(0),  
70   fHistScalevsCent(0),  
71   fHistDeltaScalevsCent(0), 
72   fHistScaleEmcalvsCent(0),      
73   fHistScale2EmcalvsCent(0),
74   fHistDeltaScale2EmcalvsCent(0),  
75   fHistChScalevsCent(0),          
76   fHistChScale2EmcalvsCent(0),   
77   fHistPtTPCvsNtrack(0), 
78   fHistPtEMCALvsNtrack(0), 
79   fHistEtvsNtrack(0),  
80   fHistScalevsNtrack(0),  
81   fHistDeltaScalevsNtrack(0),
82   fHistScaleEmcalvsNtrack(0),      
83   fHistScale2EmcalvsNtrack(0),     
84   fHistChScalevsNtrack(0),          
85   fHistChScale2EmcalvsNtrack(0),   
86   fHistTrackPtvsCent(0), 
87   fHistClusterPtvsCent(0),
88   fHistTrackEtaPhi(0), 
89   fHistClusterEtaPhi(0),
90   fHistScalevsScale2Emcal(0),      
91   fHistScalevsScaleEmcal(0),       
92   fHistScaleEmcalvsScale2Emcal(0)
93 {
94   // Constructor.
95
96   SetMakeGeneralHistograms(kTRUE);
97 }
98
99 //________________________________________________________________________
100 void AliAnalysisTaskScale::UserCreateOutputObjects()
101 {
102   // Create my user objects.
103
104   AliAnalysisTaskEmcal::UserCreateOutputObjects();
105
106   fHistPtTPCvsCent = new TH2F("fHistPtTPCvsCent", "fHistPtTPCvsCent", 101, -1, 100, 750, 0, 1500);
107   fHistPtTPCvsCent->GetXaxis()->SetTitle("Centrality (%)");
108   fHistPtTPCvsCent->GetYaxis()->SetTitle("#sum p_{T,track}^{TPC} GeV/c");
109   fHistPtTPCvsCent->GetZaxis()->SetTitle("counts");
110   fOutput->Add(fHistPtTPCvsCent);
111
112   fHistPtEMCALvsCent = new TH2F("fHistPtEMCALvsCent", "fHistPtEMCALvsCent", 101, -1, 100, 250, 0, 500);
113   fHistPtEMCALvsCent->GetXaxis()->SetTitle("Centrality (%)");
114   fHistPtEMCALvsCent->GetYaxis()->SetTitle("#sum p_{T,track}^{EMCal} GeV/c");
115   fHistPtEMCALvsCent->GetZaxis()->SetTitle("counts");
116   fOutput->Add(fHistPtEMCALvsCent);
117
118   fHistEtvsCent = new TH2F("fHistEtvsCent", "fHistEtvsCent", 101, -1, 100, 250, 0, 500);
119   fHistEtvsCent->GetXaxis()->SetTitle("Centrality (%)");
120   fHistEtvsCent->GetYaxis()->SetTitle("#sum E_{T,cluster} GeV");
121   fHistEtvsCent->GetZaxis()->SetTitle("counts");
122   fOutput->Add(fHistEtvsCent);
123
124   fHistScalevsCent = new TH2F("fHistScalevsCent", "fHistScalevsCent", 101, -1, 100, 500, 0, 5);
125   fHistScalevsCent->GetXaxis()->SetTitle("Centrality (%)");
126   fHistScalevsCent->GetYaxis()->SetTitle("s_{TPC} = (#sum E_{T,cluster} + #sum p_{T,track}^{TPC}) / #sum p_{T,track}^{TPC}");
127   fHistScalevsCent->GetZaxis()->SetTitle("counts");
128   fOutput->Add(fHistScalevsCent);
129
130   fHistDeltaScalevsCent = new TH2F("fHistDeltaScalevsCent", "fHistDeltaScalevsCent", 101, -1, 100, 500, -2.5, 2.5);
131   fHistDeltaScalevsCent->GetXaxis()->SetTitle("Centrality (%)");
132   fHistDeltaScalevsCent->GetYaxis()->SetTitle("s_{TPC}-s^{old}");
133   fHistDeltaScalevsCent->GetZaxis()->SetTitle("counts");
134   fOutput->Add(fHistDeltaScalevsCent);
135
136   fHistScaleEmcalvsCent= new TH2F("fHistScaleEmcalvsCent", "fHistScaleEmcalvsCent", 101, -1, 100, 500, 0, 5);
137   fHistScaleEmcalvsCent->GetXaxis()->SetTitle("Centrality (%)");
138   fHistScaleEmcalvsCent->GetYaxis()->SetTitle("s_{EMC}");
139   fHistScaleEmcalvsCent->GetZaxis()->SetTitle("counts");
140   fOutput->Add(fHistScaleEmcalvsCent);
141
142   fHistScale2EmcalvsCent = new TH2F("fHistScale2EmcalvsCent", "fHistScale2EmcalvsCent", 101, -1, 100, 500, 0, 5);
143   fHistScale2EmcalvsCent->GetXaxis()->SetTitle("Centrality (%)");
144   fHistScale2EmcalvsCent->GetYaxis()->SetTitle("s_{2 #times EMC}");
145   fHistScale2EmcalvsCent->GetZaxis()->SetTitle("counts");
146   fOutput->Add(fHistScale2EmcalvsCent);
147
148   fHistDeltaScale2EmcalvsCent = new TH2F("fHistDeltaScale2EmcalvsCent", "fHistDeltaScale2EmcalvsCent", 101, -1, 100, 500, -2.5, 2.5);
149   fHistDeltaScale2EmcalvsCent->GetXaxis()->SetTitle("Centrality (%)");
150   fHistDeltaScale2EmcalvsCent->GetYaxis()->SetTitle("s_{2 #times EMC}-s^{old}");
151   fHistDeltaScale2EmcalvsCent->GetZaxis()->SetTitle("counts");
152   fOutput->Add(fHistDeltaScale2EmcalvsCent);
153
154   fHistChScalevsCent = new TH2F("fHistChScalevsCent", "fHistChScalevsCent", 101, -1, 100, 500, 0, 5);
155   fHistChScalevsCent->GetXaxis()->SetTitle("Centrality (%)");
156   fHistChScalevsCent->GetYaxis()->SetTitle("s_{TPC}^{ch}");
157   fHistChScalevsCent->GetZaxis()->SetTitle("counts");
158   fOutput->Add(fHistChScalevsCent);
159
160   fHistChScale2EmcalvsCent = new TH2F("fHistChScale2EmcalvsCent", "fHistChScale2EmcalvsCent", 101, -1, 100, 500, 0, 5);
161   fHistChScale2EmcalvsCent->GetXaxis()->SetTitle("Centrality (%)");
162   fHistChScale2EmcalvsCent->GetYaxis()->SetTitle("s_{2 #times EMC}^{ch}");
163   fHistChScale2EmcalvsCent->GetZaxis()->SetTitle("counts");
164   fOutput->Add(fHistChScale2EmcalvsCent);
165
166   fHistPtTPCvsNtrack = new TH2F("fHistPtTPCvsNtrack", "fHistPtTPCvsNtrack", 800, 0, 4000,  750, 0, 1500);
167   fHistPtTPCvsNtrack->GetXaxis()->SetTitle("No. of tracks");
168   fHistPtTPCvsNtrack->GetYaxis()->SetTitle("#sum p_{T,track}^{TPC}");
169   fHistPtTPCvsNtrack->GetZaxis()->SetTitle("counts");
170   fOutput->Add(fHistPtTPCvsNtrack);
171
172   fHistPtEMCALvsNtrack = new TH2F("fHistPtEMCALvsNtrack", "fHistPtEMCALvsNtrack", 800, 0, 4000,  500, 0, 1000);
173   fHistPtEMCALvsNtrack->GetXaxis()->SetTitle("No. of tracks");
174   fHistPtEMCALvsNtrack->GetYaxis()->SetTitle("#sum p_{T,track}^{EMCal}");
175   fHistPtEMCALvsNtrack->GetZaxis()->SetTitle("counts");
176   fOutput->Add(fHistPtEMCALvsNtrack);
177
178   fHistEtvsNtrack = new TH2F("fHistEtvsNtrack", "fHistEtvsNtrack", 800,  0, 4000, 500, 0, 1000);
179   fHistEtvsNtrack->GetXaxis()->SetTitle("No. of tracks");
180   fHistEtvsNtrack->GetYaxis()->SetTitle("#sum E_{T,cluster}");
181   fHistEtvsNtrack->GetZaxis()->SetTitle("counts");
182   fOutput->Add(fHistEtvsNtrack);
183
184   fHistScalevsNtrack = new TH2F("fHistScalevsNtrack", "fHistScalevsNtrack", 800, 0, 4000,  500, 0, 5);
185   fHistScalevsNtrack->GetXaxis()->SetTitle("No. of tracks");
186   fHistScalevsNtrack->GetYaxis()->SetTitle("s_{TPC}");
187   fHistScalevsNtrack->GetZaxis()->SetTitle("counts");
188   fOutput->Add(fHistScalevsNtrack);
189
190   fHistDeltaScalevsNtrack = new TH2F("fHistDeltaScalevsNtrack", "fHistDeltaScalevsNtrack", 800, 0, 4000, 500, -2.5, 2.5);
191   fHistDeltaScalevsNtrack->GetXaxis()->SetTitle("No. of tracks");
192   fHistDeltaScalevsNtrack->GetYaxis()->SetTitle("s_{TPC}-s^{old}");
193   fHistDeltaScalevsNtrack->GetZaxis()->SetTitle("counts");
194   fOutput->Add(fHistDeltaScalevsNtrack);
195
196   fHistScaleEmcalvsNtrack = new TH2F("fHistScaleEmcalvsNtrack", "fHistScaleEmcalvsNtrack", 800, 0, 4000, 500, 0, 5);
197   fHistScaleEmcalvsNtrack->GetXaxis()->SetTitle("No. of tracks");
198   fHistScaleEmcalvsNtrack->GetYaxis()->SetTitle("s_{EMC}");
199   fHistScaleEmcalvsNtrack->GetZaxis()->SetTitle("counts");
200   fOutput->Add(fHistScaleEmcalvsNtrack);
201
202   fHistScale2EmcalvsNtrack = new TH2F("fHistScale2EmcalvsNtrack","fHistScale2EmcalvsNtrack", 800, 0, 4000, 500, 0, 5);
203   fHistScale2EmcalvsNtrack->GetXaxis()->SetTitle("No. of tracks");
204   fHistScale2EmcalvsNtrack->GetYaxis()->SetTitle("s_{2 #times EMC}");
205   fHistScale2EmcalvsNtrack->GetZaxis()->SetTitle("counts");
206   fOutput->Add(fHistScale2EmcalvsNtrack);
207
208   fHistChScalevsNtrack = new TH2F("fHistChScalevsNtrack", "fHistChScalevsNtrack", 800, 0, 4000, 500, 0, 5);
209   fHistChScalevsNtrack->GetXaxis()->SetTitle("No. of tracks");
210   fHistChScalevsNtrack->GetYaxis()->SetTitle("s_{TPC}^{ch}");
211   fHistChScalevsNtrack->GetZaxis()->SetTitle("counts");
212   fOutput->Add(fHistChScalevsNtrack);
213
214   fHistChScale2EmcalvsNtrack = new TH2F("fHistChScale2EmcalvsNtrack", "fHistChScale2EmcalvsNtrack", 800,  0, 4000, 500, 0, 5);
215   fHistChScale2EmcalvsNtrack->GetXaxis()->SetTitle("No. of tracks");
216   fHistChScale2EmcalvsNtrack->GetYaxis()->SetTitle("s_{2 #times EMC}^{ch}");
217   fHistChScale2EmcalvsNtrack->GetZaxis()->SetTitle("counts");
218   fOutput->Add(fHistChScale2EmcalvsNtrack);
219
220   fHistTrackPtvsCent = new TH2F("fHistTrackPtvsCent", "fHistTrackPtvsCent", 101, -1, 100, 500, 0, 100);
221   fHistTrackPtvsCent->GetXaxis()->SetTitle("Centrality (%)");
222   fHistTrackPtvsCent->GetYaxis()->SetTitle("p_{T,track} GeV/c");
223   fHistTrackPtvsCent->GetZaxis()->SetTitle("counts");
224   fOutput->Add(fHistTrackPtvsCent);
225
226   fHistClusterPtvsCent = new TH2F("fHistClusterPtvsCent", "fHistClusterPtvsCent", 101, -1, 100, 500, 0, 100);
227   fHistClusterPtvsCent->GetXaxis()->SetTitle("Centrality (%)");
228   fHistClusterPtvsCent->GetYaxis()->SetTitle("E_{T,cluster} GeV");
229   fHistClusterPtvsCent->GetZaxis()->SetTitle("counts");
230   fOutput->Add(fHistClusterPtvsCent);
231
232   fHistTrackEtaPhi = new TH2F("fHistTrackEtaPhi", "fHistTrackEtaPhi", 100, -1.0, 1.0, 101, 0, 2.02*TMath::Pi());
233   fHistTrackEtaPhi->GetXaxis()->SetTitle("#eta");
234   fHistTrackEtaPhi->GetYaxis()->SetTitle("#phi");
235   fHistTrackEtaPhi->GetZaxis()->SetTitle("counts");
236   fOutput->Add(fHistTrackEtaPhi);
237
238   fHistClusterEtaPhi = new TH2F("fHistClusterEtaPhi", "fHistClusterEtaPhi", 100, -1.0, 1.0, 101, 0, 2.02*TMath::Pi());
239   fHistClusterEtaPhi->GetXaxis()->SetTitle("#eta");
240   fHistClusterEtaPhi->GetYaxis()->SetTitle("#phi");
241   fHistClusterEtaPhi->GetZaxis()->SetTitle("counts");
242   fOutput->Add(fHistClusterEtaPhi);
243
244   fHistScalevsScale2Emcal = new TH2F("fHistScalevsScale2Emcal", "fHistScalevsScale2Emcal",500, 0, 5, 500, 0, 5);
245   fHistScalevsScale2Emcal->GetXaxis()->SetTitle("s_{TPC}");
246   fHistScalevsScale2Emcal->GetYaxis()->SetTitle("s_{2 #times EMC}");
247   fHistScalevsScale2Emcal->GetZaxis()->SetTitle("counts");
248   fOutput->Add(fHistScalevsScale2Emcal);
249
250   fHistScalevsScaleEmcal = new TH2F("fHistScalevsScaleEmcal", "fHistScalevsScaleEmcal", 500, 0, 5, 500, 0, 5);
251   fHistScalevsScaleEmcal->GetXaxis()->SetTitle("s_{TPC}");
252   fHistScalevsScaleEmcal->GetYaxis()->SetTitle("s_{EMC}");
253   fHistScalevsScaleEmcal->GetZaxis()->SetTitle("counts");
254   fOutput->Add(fHistScalevsScaleEmcal);
255
256   fHistScaleEmcalvsScale2Emcal = new TH2F("fHistScaleEmcalvsScale2Emcal", "fHistScaleEmcalvsScale2Emcal", 500, 0, 5, 500, 0, 5);
257   fHistScaleEmcalvsScale2Emcal->GetXaxis()->SetTitle("s_{EMC}");
258   fHistScaleEmcalvsScale2Emcal->GetYaxis()->SetTitle("s_{2 #times EMC}");
259   fHistScaleEmcalvsScale2Emcal->GetZaxis()->SetTitle("counts");
260   fOutput->Add(fHistScaleEmcalvsScale2Emcal);
261
262   PostData(1, fOutput);
263 }
264
265 //________________________________________________________________________
266 Double_t AliAnalysisTaskScale::GetScaleFactor(Double_t cent)
267 {
268   // Get scale function.
269
270   Double_t scale = -1;
271   if (fScaleFunction)
272     scale = fScaleFunction->Eval(cent);
273   return scale;
274 }
275
276 //________________________________________________________________________
277 Bool_t AliAnalysisTaskScale::FillHistograms() 
278 {
279   // Execute on each event.
280
281   const Double_t EmcalMinPhi = fGeom->GetArm1PhiMin() * TMath::DegToRad();
282   const Double_t EmcalMaxPhi = fGeom->GetArm1PhiMax() * TMath::DegToRad();
283   const Double_t EmcalWidth = (EmcalMaxPhi-EmcalMinPhi)/2.0;
284
285   Double_t ptTPC   = 0;
286   Double_t ptEMCAL = 0;
287   Double_t ptEMCAL2 = 0;
288
289   const Int_t Ntracks = fTracks->GetEntries();
290   for (Int_t iTracks = 0; iTracks < Ntracks; ++iTracks) {
291     AliVParticle *track = static_cast<AliVParticle*>(fTracks->At(iTracks));
292
293     if (!track)
294       continue;
295
296     if (!AcceptTrack(track))
297       continue;
298
299     if (TMath::Abs(track->Eta()) > 0.7)   // only accept tracks in the EMCal eta range
300       continue;
301
302     fHistTrackPtvsCent->Fill(fCent,track->Pt());
303     fHistTrackEtaPhi->Fill(track->Eta(),track->Phi());
304     ptTPC += track->Pt();
305     if ((track->Phi() > (EmcalMaxPhi+EmcalWidth)) || (track->Phi() < (EmcalMinPhi-EmcalWidth))) continue;
306     ptEMCAL2 += track->Pt();
307     if ((track->Phi() > EmcalMaxPhi) || (track->Phi() < EmcalMinPhi)) continue;
308     ptEMCAL += track->Pt();
309   }
310
311   if (ptTPC == 0) 
312     return kFALSE;
313   
314   Double_t Et = 0;
315   const Int_t Nclus = fCaloClusters->GetEntries();
316   for (Int_t iClus = 0; iClus < Nclus; ++iClus) {
317     AliVCluster *c = static_cast<AliVCluster*>(fCaloClusters->At(iClus));
318     if (!c)
319       continue;
320
321     if (!AcceptCluster(c))
322       continue;
323
324     TLorentzVector nPart;
325     c->GetMomentum(nPart, fVertex);
326
327     fHistClusterPtvsCent->Fill(fCent, nPart.Pt());
328     fHistClusterEtaPhi->Fill(nPart.Eta(), nPart.Phi());
329
330     Et += nPart.Pt();
331   }
332  
333   Double_t scalecalc         = -1;
334   if (ptEMCAL > 0 && Et > 0 && ptTPC > 0)
335     scalecalc         =  ((Et + ptEMCAL) / fEmcalArea) * (fTpcArea / ptTPC);
336   const Double_t scale             = GetScaleFactor(fCent);
337   Double_t scalecalcemcal          = -1;
338   if (ptEMCAL > 0)
339     scalecalcemcal                 = (Et+ptEMCAL)/ptEMCAL;
340   Double_t scalecalcemcal2         = -1;
341   Double_t Chscalecalcemcal2       = -1;
342   if (ptEMCAL2 > 0){
343     scalecalcemcal2                = 2*(Et+ptEMCAL)/ptEMCAL2;
344     Chscalecalcemcal2              = 2*ptEMCAL/ptEMCAL2;}
345   const Double_t Chscalecalcemcal  = ((ptEMCAL) / fEmcalArea) * (fTpcArea / ptTPC);
346
347   fHistScaleEmcalvsCent->Fill(fCent,scalecalcemcal);      
348   fHistScale2EmcalvsCent->Fill(fCent,scalecalcemcal2);     
349   fHistChScalevsCent->Fill(fCent,Chscalecalcemcal);    
350   fHistChScale2EmcalvsCent->Fill(fCent,Chscalecalcemcal2);   
351   fHistScaleEmcalvsNtrack->Fill(Ntracks,scalecalcemcal);      
352   fHistScale2EmcalvsNtrack->Fill(Ntracks,scalecalcemcal2);     
353   fHistChScalevsNtrack->Fill(Ntracks,Chscalecalcemcal);    
354   fHistChScale2EmcalvsNtrack->Fill(Ntracks,Chscalecalcemcal2);   
355   fHistPtTPCvsCent->Fill(fCent, ptTPC);
356   fHistPtEMCALvsCent->Fill(fCent, ptEMCAL);
357   fHistEtvsCent->Fill(fCent, Et);
358   fHistScalevsCent->Fill(fCent, scalecalc);
359   fHistDeltaScalevsCent->Fill(fCent, scalecalc - scale);
360   fHistDeltaScale2EmcalvsCent->Fill(fCent, scalecalcemcal2 - scale);
361   fHistPtTPCvsNtrack->Fill(Ntracks, ptTPC);
362   fHistPtEMCALvsNtrack->Fill(Ntracks, ptEMCAL);
363   fHistEtvsNtrack->Fill(Ntracks, Et);
364   fHistScalevsNtrack->Fill(Ntracks, scalecalc);
365   fHistDeltaScalevsNtrack->Fill(Ntracks, scalecalc - scale);
366   fHistScalevsScale2Emcal->Fill(scalecalc,scalecalcemcal2);      
367   fHistScalevsScaleEmcal->Fill(scalecalc,scalecalcemcal); 
368   fHistScaleEmcalvsScale2Emcal->Fill(scalecalcemcal,scalecalcemcal2);
369
370   return kTRUE;
371 }
372
373
374 //________________________________________________________________________
375 void AliAnalysisTaskScale::ExecOnce() 
376 {
377   AliAnalysisTaskEmcal::ExecOnce();
378
379   const Double_t EmcalMinEta = fGeom->GetArm1EtaMin();
380   const Double_t EmcalMaxEta = fGeom->GetArm1EtaMax();
381   const Double_t EmcalMinPhi = fGeom->GetArm1PhiMin() * TMath::DegToRad();
382   const Double_t EmcalMaxPhi = fGeom->GetArm1PhiMax() * TMath::DegToRad();
383
384   fEmcalArea  = (EmcalMaxPhi - EmcalMinPhi) * (EmcalMinEta - EmcalMaxEta);
385
386   AliParticleContainer *partCont = GetParticleContainer(0);
387   if (!partCont) {
388     AliError(Form("%s: No particle container found! Assuming tpc area = 1...",GetName()));
389     fTpcArea = 1;
390     return;
391   }
392
393   Float_t TpcMaxPhi = partCont->GetParticlePhiMax();
394   Float_t TpcMinPhi = partCont->GetParticlePhiMin();
395   
396   if (TpcMaxPhi > TMath::Pi()*2) TpcMaxPhi = TMath::Pi()*2;
397   if (TpcMinPhi < 0) TpcMinPhi = 0;
398
399   fTpcArea = (TpcMaxPhi - TpcMinPhi) * (EmcalMinEta - EmcalMaxEta);
400 }