]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGGA/EMCALTasks/AliAnalysisTaskSATR.cxx
including setter to enable debugging printouts
[u/mrichter/AliRoot.git] / PWGGA / EMCALTasks / AliAnalysisTaskSATR.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 // --- Root ---
19 #include <TChain.h>
20 #include <TFile.h>
21 #include <TH1F.h>
22 #include <TH2F.h>
23 #include <TList.h>
24 #include <TClonesArray.h>
25
26 // --- AliRoot ---
27 #include "AliAnalysisTaskSE.h"
28 #include "AliAnalysisManager.h"
29 #include "AliVEvent.h"
30 #include "AliVCaloTrigger.h"
31 #include "AliVCluster.h"
32 #include "AliEMCALGeometry.h"
33 #include "AliLog.h"
34 #include "AliCaloCalibPedestal.h"
35 #include "AliCDBEntry.h"
36 #include "AliCDBManager.h"
37 #include "AliAnalysisTaskEMCALClusterizeFast.h"
38 #include "AliCentrality.h"
39
40 #include "AliAnalysisTaskSATR.h"
41
42 ClassImp(AliAnalysisTaskSATR)
43
44 //________________________________________________________________________
45 AliAnalysisTaskSATR::AliAnalysisTaskSATR() : 
46   AliAnalysisTaskSE("AliAnalysisTaskSATR"),
47   fL0Calib(0.065),
48   fCaloClustersName("CaloClusters"),
49   fTriggerClustersName("triggerClusters"),
50   fMinCutL0Amp(-1),
51   fMaxCutL0Amp(999),
52   fMinCutClusEnergy(-1),
53   fMaxCutClusEnergy(999),
54   fTimeCutOn(0),
55   fMinL0Time(-20),
56   fMaxL0Time(20),
57   fMinClusTime(-1),
58   fMaxClusTime(1),
59   fCheckDeadClusters(0),
60   fPedestal(0),
61   fLoadPed(0),
62   fOCDBpath(""),
63   fMinDistanceFromBadTower(0.5),
64   fClusterizer(0),
65   fTriggerClusterizer(0),
66   fGeom(0),
67   fRun(-1),
68   fOutput(0),
69   fHistEclus(0),
70   fHistEmaxClus(0),
71   fHistEtavsPhiMaxClus(0),
72   fHistEtavsEmaxClus(0),
73   fHistPhivsEmaxClus(0),
74   fHistTOFvsEclus(0),
75   fHistTOFvsEclusC(0),
76   fHistNcellsvsEclus(0),
77   fHistAmpTClus(0),
78   fHistAmpMaxTClus(0),
79   fHistEtavsPhiMaxTClus(0),
80   fHistEmaxClusvsAmpMaxTClus(0),
81   fHistEmaxClusvsAmpMatchedTClus(0),
82   fHistEmaxClusNotMatchingTClus(0),
83   fHistEtavsPhiMaxClusNotMatchingTClus(0),
84   fHistEmatchedClusvsAmpMaxTClus(0),
85   fHistAmpMaxTClusNotMatchingClus(0),
86   fHistEtavsPhiMaxTClusNotMatchingClus(0),
87   fHistIdxMaxClusvsIdxMaxTClus(0),
88   fHistPhiMaxClusvsPhiMaxTClus(0),
89   fHistEtaMaxClusvsEtaMaxTClus(0),
90   fHistTOFmaxClusvsTimeMaxTClus(0),
91   fHistEmatchedClusvsAmpMatchedTClus(0),
92   fHistEmatchedClus(0),
93   fHistEmaxMatchedClus(0),
94   fHistAmpL1TimeSum(0),
95   fHistAmpMaxL1TimeSum(0),
96   fHistAmpMaxL1TimeSumVScent(0),
97   fHistAmpFastORvsAmpL1TimeSum(0),
98   fHistAmpFastOR(0),
99   fHistAmpMaxFastOR(0),
100   fHistTimeFastOR(0),
101   fHistEtavsPhiFastOR(0),
102   fHistEtavsPhiMaxFastOR(0),
103   fHistTimeDispFastOR(0),
104   fHistTimevsL0TimeFastOR(0),
105   fHistNtimesFastOR(0),
106   fHistEcells(0),
107   fHistEmaxCell(0),
108   fHistTOFvsEcells(0),
109   fHistTOFvsEcellsC(0),
110   fHistEmaxCellvsAmpFastOR(0)
111 {
112   // Constructor
113 }
114
115 //________________________________________________________________________
116 AliAnalysisTaskSATR::AliAnalysisTaskSATR(const char *name) : 
117   AliAnalysisTaskSE(name),
118   fL0Calib(0.065),
119   fCaloClustersName("CaloClusters"),
120   fTriggerClustersName("triggerClusters"),
121   fMinCutL0Amp(-1),
122   fMaxCutL0Amp(999),
123   fMinCutClusEnergy(-1),
124   fMaxCutClusEnergy(999),
125   fTimeCutOn(0),
126   fMinL0Time(-20),
127   fMaxL0Time(20),
128   fMinClusTime(-1),
129   fMaxClusTime(1),
130   fCheckDeadClusters(0),
131   fPedestal(0),
132   fLoadPed(0),
133   fOCDBpath(""),
134   fMinDistanceFromBadTower(0.5),
135   fClusterizer(0),
136   fTriggerClusterizer(0),
137   fGeom(0),
138   fRun(-1),
139   fOutput(0),
140   fHistEclus(0),
141   fHistEmaxClus(0),
142   fHistEtavsPhiMaxClus(0),
143   fHistEtavsEmaxClus(0),
144   fHistPhivsEmaxClus(0),
145   fHistTOFvsEclus(0),
146   fHistTOFvsEclusC(0),
147   fHistNcellsvsEclus(0),
148   fHistAmpTClus(0),
149   fHistAmpMaxTClus(0),
150   fHistEtavsPhiMaxTClus(0),
151   fHistEmaxClusvsAmpMaxTClus(0),
152   fHistEmaxClusvsAmpMatchedTClus(0),
153   fHistEmaxClusNotMatchingTClus(0),
154   fHistEtavsPhiMaxClusNotMatchingTClus(0),
155   fHistEmatchedClusvsAmpMaxTClus(0),
156   fHistAmpMaxTClusNotMatchingClus(0),
157   fHistEtavsPhiMaxTClusNotMatchingClus(0),
158   fHistIdxMaxClusvsIdxMaxTClus(0),
159   fHistPhiMaxClusvsPhiMaxTClus(0),
160   fHistEtaMaxClusvsEtaMaxTClus(0),
161   fHistTOFmaxClusvsTimeMaxTClus(0),
162   fHistEmatchedClusvsAmpMatchedTClus(0),
163   fHistEmatchedClus(0),
164   fHistEmaxMatchedClus(0),
165   fHistAmpL1TimeSum(0),
166   fHistAmpMaxL1TimeSum(0),
167   fHistAmpMaxL1TimeSumVScent(0),
168   fHistAmpFastORvsAmpL1TimeSum(0),
169   fHistAmpFastOR(0),
170   fHistAmpMaxFastOR(0),
171   fHistTimeFastOR(0),
172   fHistEtavsPhiFastOR(0),
173   fHistEtavsPhiMaxFastOR(0),
174   fHistTimeDispFastOR(0),
175   fHistTimevsL0TimeFastOR(0),
176   fHistNtimesFastOR(0),
177   fHistEcells(0),
178   fHistEmaxCell(0),
179   fHistTOFvsEcells(0),
180   fHistTOFvsEcellsC(0),
181   fHistEmaxCellvsAmpFastOR(0)
182 {
183   // Constructor
184
185   DefineOutput(1, TList::Class()); 
186 }
187
188 //________________________________________________________________________
189 AliAnalysisTaskSATR::~AliAnalysisTaskSATR()
190 {
191
192 }
193
194 //________________________________________________________________________
195 void AliAnalysisTaskSATR::UserCreateOutputObjects()
196 {
197   // Create histograms
198   
199   const Int_t     L0Ampbins   = 100;
200   const Float_t   L0Amplow    = 0;
201   const Float_t   L0Ampup     = 100;
202   const Int_t     L0Timebins  = 20;
203   const Float_t   L0Timelow   = 0;
204   const Float_t   L0Timeup    = 20;
205   const Int_t     Ebins       = 100;
206   const Float_t   Elow        = 0;
207   const Float_t   Eup         = 25;
208   const Int_t     TOFbins     = 100;
209   const Float_t   TOFlow      = 0;
210   const Float_t   TOFup       = 1e-6;
211   const Float_t   L1Amplow    = 0;
212   const Float_t   L1Ampup     = 400;
213   const Int_t     L1Ampbins   = 400;
214   const Int_t     Indexesbins = 1440;
215   const Int_t     Indexeslow  = 0;
216   const Int_t     Indexesup   = 2880;
217   const Int_t     nPhibins    = 60;
218   const Int_t     nPhilow     = 0;
219   const Int_t     nPhiup      = 60;
220   const Int_t     nEtabins    = 48;
221   const Int_t     nEtalow     = 0;
222   const Int_t     nEtaup      = 48;
223   const Int_t     RowTrgbins  = 60;
224   const Int_t     RowTrglow   = 0;
225   const Int_t     RowTrgup    = 60;
226   const Int_t     ColTrgbins  = 48;
227   const Int_t     ColTrglow   = 0;
228   const Int_t     ColTrgup    = 48;
229   
230   if (fClusterizer)
231     fCaloClustersName = fClusterizer->GetNewClusterArrayName();
232   
233   if (fTriggerClusterizer)
234     fTriggerClustersName = fTriggerClusterizer->GetNewClusterArrayName();
235   
236   fOutput = new TList();
237   fOutput->SetOwner();
238
239   fHistEclus = new TH1F("fHistEclus","Energy Cluster spectrum", Ebins, Elow, Eup);
240   fHistEclus->GetXaxis()->SetTitle("Cluster Energy [GeV]");
241   fHistEclus->GetYaxis()->SetTitle("counts");
242   
243   fHistEmaxClus = new TH1F("fHistEmaxClus","Maximum Cluster Energy per event spectrum",Ebins, Elow, Eup);
244   fHistEmaxClus->GetXaxis()->SetTitle("Maximum Cluster Energy [GeV]");
245   fHistEmaxClus->GetYaxis()->SetTitle("counts");
246   
247   fHistEtavsPhiMaxClus = new TH2I("fHistEtavsPhiMaxClus","nEta vs nPhi of maximum cluster per event", nEtabins, nEtalow, nEtaup, nPhibins, nPhilow, nPhiup);
248   fHistEtavsPhiMaxClus->GetXaxis()->SetTitle("nEta");
249   fHistEtavsPhiMaxClus->GetYaxis()->SetTitle("nPhi");
250         
251   fHistEtavsEmaxClus = new TH2F("fHistEtavsEmaxClus","nEta vs Energy of maximum cluster per event", nEtabins, nEtalow, nEtaup, Ebins, Elow, Eup);
252   fHistEtavsEmaxClus->GetYaxis()->SetTitle("Energy [GeV]");
253   fHistEtavsEmaxClus->GetXaxis()->SetTitle("nEta");
254   
255   fHistPhivsEmaxClus = new TH2F("fHistPhivsEmaxClus","nPhi vs Energy of maximum cluster per event", nPhibins, nPhilow, nPhiup, Ebins, Elow, Eup);
256   fHistPhivsEmaxClus->GetYaxis()->SetTitle("Energy [GeV]");
257   fHistPhivsEmaxClus->GetXaxis()->SetTitle("nPhi");
258   
259   fHistTOFvsEclus = new TH2F("fHistTOFvsEclus","TOF vs Energy of clusters", TOFbins, TOFlow, TOFup, Ebins, Elow, Eup);
260   fHistTOFvsEclus->GetXaxis()->SetTitle("TOF [s]");
261   fHistTOFvsEclus->GetYaxis()->SetTitle("Energy [GeV]");
262         
263   fHistTOFvsEclusC = new TH2F("fHistTOFvsEclusC","TOF vs Energy of clusters (corrected)", TOFbins, TOFlow, TOFup, Ebins, Elow, Eup);
264   fHistTOFvsEclusC->GetXaxis()->SetTitle("TOF");
265   fHistTOFvsEclusC->GetYaxis()->SetTitle("E [GeV]");
266   
267   fHistNcellsvsEclus = new TH2F("fHistNcellsvsEclus","N_{cells} vs Energy of clusters", 20, 0, 20, Ebins, Elow, Eup);
268   fHistNcellsvsEclus->GetXaxis()->SetTitle("N_{cells}");
269   fHistNcellsvsEclus->GetYaxis()->SetTitle("E [GeV]");
270   
271   fHistAmpTClus = new TH1F("fHistAmpTClus","L0 Trigger Candidate Amplitude spectrum", L0Ampbins, L0Amplow, L0Ampup * 4);
272   fHistAmpTClus->GetXaxis()->SetTitle("L0 Trigger Candidate Amplitude");
273   fHistAmpTClus->GetYaxis()->SetTitle("counts");
274         
275   fHistAmpMaxTClus = new TH1F("fHistAmpMaxTClus","Maximum L0 Trigger Candidate Amplitude per event spectrum", L0Ampbins, L0Amplow, L0Ampup * 4);
276   fHistAmpMaxTClus->GetXaxis()->SetTitle("Maximum L0 Trigger Candidate Amplitude");
277   fHistAmpMaxTClus->GetYaxis()->SetTitle("counts");
278   
279   fHistEtavsPhiMaxTClus = new TH2I("fHistEtavsPhiMaxTClus","nEta vs nPhi of maximum L0 Trigger Candidate per event", nEtabins, nEtalow, nEtaup, nPhibins, nPhilow, nPhiup);
280   fHistEtavsPhiMaxTClus->GetXaxis()->SetTitle("nEta");
281   fHistEtavsPhiMaxTClus->GetYaxis()->SetTitle("nPhi");
282   
283   fHistEmaxClusvsAmpMaxTClus = new TH2F("fHistEmaxClusvsAmpMaxTClus","Maximum Cluster Energy vs Maximum L0 Trigger Candidate Amplitude", L0Ampbins, L0Amplow, L0Ampup * 4, Ebins, Elow, Eup);
284   fHistEmaxClusvsAmpMaxTClus->GetXaxis()->SetTitle("Maximum L0 Trigger Candidate Amplitude");
285   fHistEmaxClusvsAmpMaxTClus->GetYaxis()->SetTitle("Maximum Cluster Energy [GeV]");
286   
287   fHistEmaxClusvsAmpMatchedTClus = new TH2F("fHistEmaxClusvsAmpMatchedTClus","Maximum Cluster Energy vs Matched L0 Trigger Candidate Amplitude", L0Ampbins, L0Amplow, L0Ampup * 4, Ebins, Elow, Eup);
288   fHistEmaxClusvsAmpMatchedTClus->GetXaxis()->SetTitle("Matched L0 Trigger Candidate Amplitude");
289   fHistEmaxClusvsAmpMatchedTClus->GetYaxis()->SetTitle("Maximum Cluster Energy [GeV]");
290
291   fHistEmaxClusNotMatchingTClus = new TH1F("fHistEmaxClusNotMatchingTClus","Maximum Cluster Energy per event not matching any fake ALTRO data", Ebins, Elow, Eup);
292   fHistEmaxClusNotMatchingTClus->GetXaxis()->SetTitle("Maximum Cluster Energy [GeV]");
293   fHistEmaxClusNotMatchingTClus->GetYaxis()->SetTitle("counts");
294
295   fHistEtavsPhiMaxClusNotMatchingTClus = new TH2I("fHistEtavsPhiMaxClusNotMatchingTClus","nEta vs nPhi of maximum cluster per event not matching any fake ALTRO data", nEtabins, nEtalow, nEtaup, nPhibins, nPhilow, nPhiup);
296   fHistEtavsPhiMaxClusNotMatchingTClus->GetXaxis()->SetTitle("nEta");
297   fHistEtavsPhiMaxClusNotMatchingTClus->GetYaxis()->SetTitle("nPhi");
298   
299   fHistEmatchedClusvsAmpMaxTClus = new TH2F("fHistEmatchedClusvsAmpMaxTClus","Matched Cluster Energy vs Maximum L0 Trigger Candidate Amplitude", L0Ampbins, L0Amplow, L0Ampup * 4, Ebins, Elow, Eup);
300   fHistEmatchedClusvsAmpMaxTClus->GetXaxis()->SetTitle("Maximum L0 Trigger Candidate Amplitude");
301   fHistEmatchedClusvsAmpMaxTClus->GetYaxis()->SetTitle("Matched Cluster Energy [GeV]");
302   
303   fHistAmpMaxTClusNotMatchingClus = new TH1F("fHistAmpMaxTClusNotMatchingClus","Maximum L0 Trigger Candidate Amplitude not matching any FEE data", L0Ampbins, L0Amplow, L0Ampup * 4);
304   fHistAmpMaxTClusNotMatchingClus->GetXaxis()->SetTitle("Maximum L0 Trigger Candidate Amplitude");
305   fHistAmpMaxTClusNotMatchingClus->GetYaxis()->SetTitle("counts");
306   
307   fHistEtavsPhiMaxTClusNotMatchingClus = new TH2I("fHistEtavsPhiMaxTClusNotMatchingClus","nEta vs nPhi of maximum L0 Trigger Candidate Amplitude per event not matching any FEE data", nEtabins, nEtalow, nEtaup, nPhibins, nPhilow, nPhiup);
308   fHistEtavsPhiMaxTClusNotMatchingClus->GetXaxis()->SetTitle("nEta");
309   fHistEtavsPhiMaxTClusNotMatchingClus->GetYaxis()->SetTitle("nPhi");
310   
311   fHistIdxMaxClusvsIdxMaxTClus = new TH2I("fHistIdxMaxClusvsIdxMaxTClus","Maximum Cluster Index vs Maximum L0 Trigger Candidate Index", Indexesbins, Indexeslow, Indexesup, Indexesbins, Indexeslow, Indexesup);
312   fHistIdxMaxClusvsIdxMaxTClus->GetXaxis()->SetTitle("Maximum L0 Trigger Candidate Index");
313   fHistIdxMaxClusvsIdxMaxTClus->GetYaxis()->SetTitle("Maximum Cluster Index");
314   
315   fHistPhiMaxClusvsPhiMaxTClus = new TH2I("fHistPhiMaxClusvsPhiMaxTClus","Maximum Cluster nPhi vs Maximum L0 Trigger Candidate nPhi", nPhibins, nPhilow, nPhiup, nPhibins, nPhilow, nPhiup);
316   fHistPhiMaxClusvsPhiMaxTClus->GetXaxis()->SetTitle("Maximum L0 Trigger Candidate nPhi");
317   fHistPhiMaxClusvsPhiMaxTClus->GetYaxis()->SetTitle("Maximum Cluster nPhi");
318
319   fHistEtaMaxClusvsEtaMaxTClus = new TH2I("fHistEtaMaxClusvsEtaMaxTClus","Maximum Cluster nEta vs Maximum L0 Trigger Candidate nEta", nEtabins, nEtalow, nEtaup, nEtabins, nEtalow, nEtaup);
320   fHistEtaMaxClusvsEtaMaxTClus->GetXaxis()->SetTitle("Maximum L0 Trigger Candidate nEta");
321   fHistEtaMaxClusvsEtaMaxTClus->GetYaxis()->SetTitle("Maximum Cluster nEta");
322   
323   fHistTOFmaxClusvsTimeMaxTClus = new TH2F("fHistTOFmaxClusvsTimeMaxTClus","Maximum Cluster TOF vs Maximum L0 Trigger Candidate Time", L0Timebins, L0Timelow, L0Timeup, TOFbins, TOFlow, TOFup);
324   fHistTOFmaxClusvsTimeMaxTClus->GetXaxis()->SetTitle("Maximum L0 Trigger Candidate Time");
325   fHistTOFmaxClusvsTimeMaxTClus->GetYaxis()->SetTitle("Maximum Cluster TOF");
326
327   fHistEmatchedClusvsAmpMatchedTClus = new TH2F("fHistEmatchedClusvsAmpMatchedTClus","Matched Cluster Energy vs Matched L0 Trigger Candidate Amplitude", L0Ampbins, L0Amplow, L0Ampup * 4, Ebins, Elow, Eup);
328   fHistEmatchedClusvsAmpMatchedTClus->GetXaxis()->SetTitle("Matched Trigger Candidate Amplitude");
329   fHistEmatchedClusvsAmpMatchedTClus->GetYaxis()->SetTitle("Matched Cluster Energy [GeV]");
330   
331   fHistEmatchedClus = new TH1F("fHistEmatchedClus","Cluster Energy spectrum that matched L0 Trigger Candidate",Ebins, Elow, Eup);
332   fHistEmatchedClus->GetXaxis()->SetTitle("Cluster Energy [GeV]");
333   fHistEmatchedClus->GetYaxis()->SetTitle("counts");
334   
335   fHistEmaxMatchedClus = new TH1F("fHistEmaxMatchedClus","Maximum Cluster Energy per event that matched L0 Trigger Candidate",Ebins, Elow, Eup);
336   fHistEmaxMatchedClus->GetXaxis()->SetTitle("Cluster Energy [GeV]");
337   fHistEmaxMatchedClus->GetYaxis()->SetTitle("counts");
338   
339   fHistAmpL1TimeSum = new TH1F("fHistAmpL1TimeSum","L1 Time Sum Amplitude spectrum",L1Ampbins, L1Amplow, L1Ampup);
340   fHistAmpL1TimeSum->GetXaxis()->SetTitle("L1 Time Sum Amplitude");
341   fHistAmpL1TimeSum->GetYaxis()->SetTitle("counts");
342   
343   fHistAmpMaxL1TimeSum = new TH1F("fHistAmpMaxL1TimeSum","L1 Time Sum Amplitude per event spectrum",L1Ampbins, L1Amplow, L1Ampup);
344   fHistAmpMaxL1TimeSum->GetXaxis()->SetTitle("Maximum L1 Time Sum Amplitude");
345   fHistAmpMaxL1TimeSum->GetYaxis()->SetTitle("counts");
346
347   fHistAmpMaxL1TimeSumVScent = new TH2F("fHistAmpMaxL1TimeSumVScent","L1 Time Sum Amplitude vs. centrality", 100, 0, 100, L1Ampbins, L1Amplow, L1Ampup);
348   fHistAmpMaxL1TimeSumVScent->GetXaxis()->SetTitle("Centrality [%]");
349   fHistAmpMaxL1TimeSumVScent->GetYaxis()->SetTitle("Maximum L1 Time Sum Amplitude");
350
351   fHistAmpFastORvsAmpL1TimeSum = new TH2F("fHistAmpFastORvsAmpL1TimeSum","FastOR Amplitude vs L1 Time Sum Amplitude",L0Ampbins, L0Amplow, L0Ampup,L1Ampbins, L1Amplow, L1Ampup);
352   fHistAmpFastORvsAmpL1TimeSum->GetXaxis()->SetTitle("FastOR Amplitude");
353   fHistAmpFastORvsAmpL1TimeSum->GetYaxis()->SetTitle("L1 Time Sum Amplitude");
354   
355   fHistAmpFastOR = new TH1F("fHistAmpFastOR","FastOR Amplitude spectrum",L0Ampbins, L0Amplow, L0Ampup);
356   fHistAmpFastOR->GetXaxis()->SetTitle("FastOR Amplitude");
357   fHistAmpFastOR->GetYaxis()->SetTitle("counts");
358         
359   fHistAmpMaxFastOR = new TH1F("fHistAmpMaxFastOR","Maximum FastOR Amplitude per event spectrum",L0Ampbins, L0Amplow, L0Ampup);
360   fHistAmpMaxFastOR->GetXaxis()->SetTitle("max FastOR amplitude");
361   fHistAmpMaxFastOR->GetYaxis()->SetTitle("counts");
362   
363   fHistTimeFastOR = new TH1F("fHistTimeFastOR","FastOR L0 Time distribution", L0Timebins, L0Timelow, L0Timeup);
364   fHistTimeFastOR->GetXaxis()->SetTitle("FastOR L0 Time");
365   fHistTimeFastOR->GetYaxis()->SetTitle("counts");
366
367   fHistEtavsPhiFastOR = new TH2I("fHistEtavsPhiFastOR","FastOR gCol vs gRow", ColTrgbins, ColTrglow, ColTrgup,  RowTrgbins, RowTrglow, RowTrgup);
368   fHistEtavsPhiFastOR->GetXaxis()->SetTitle("FastOR gCol");
369   fHistEtavsPhiFastOR->GetYaxis()->SetTitle("FastOR gRow");
370         
371   fHistEtavsPhiMaxFastOR = new TH2I("fHistEtavsPhiMaxFastOR","Maximum FastOR gCol vs gRow", ColTrgbins, ColTrglow, ColTrgup,  RowTrgbins, RowTrglow, RowTrgup);
372   fHistEtavsPhiMaxFastOR->GetXaxis()->SetTitle("Maximum FastOR gCol");
373   fHistEtavsPhiMaxFastOR->GetYaxis()->SetTitle("Maximum FastOR gRow");
374   
375   fHistTimeDispFastOR = new TH1F("fHistTimeDispFastOR", "Dispersion from the maximum FastOR of L0 Time per event", L0Timebins * 2, -L0Timeup + .5, L0Timeup - .5);
376   fHistTimeDispFastOR->GetXaxis()->SetTitle("Disperion");
377   fHistTimeDispFastOR->GetYaxis()->SetTitle("counts");
378   
379   fHistTimevsL0TimeFastOR = new TH2F("fHistTimevsL0TimeFastOR", "FastOR Time vs FastOR L0 Time", L0Timebins, L0Timelow, L0Timeup, L0Timebins, L0Timelow, L0Timeup);
380   fHistTimevsL0TimeFastOR->GetXaxis()->SetTitle("Time");
381   fHistTimevsL0TimeFastOR->GetYaxis()->SetTitle("L0 Time");
382   
383   fHistNtimesFastOR = new TH1I("fHistNtimesFastOR", "FastOR NL0Times distribution", 5, 0, 5);
384   fHistNtimesFastOR->GetXaxis()->SetTitle("NL0Times");
385   fHistNtimesFastOR->GetYaxis()->SetTitle("counts");
386   
387   fHistEcells = new TH1F("fHistEcells","Cell Energy spectrum", Ebins, Elow, Eup);
388   fHistEcells->GetXaxis()->SetTitle("Energy [GeV]");
389   fHistEcells->GetYaxis()->SetTitle("counts");
390   
391   fHistEmaxCell = new TH1F("fHistEmaxCell","Maximum Cell Energy per event spectrum", Ebins, Elow, Eup);
392   fHistEmaxCell->GetXaxis()->SetTitle("Maximum Cell Energy [GeV]");
393   fHistEmaxCell->GetYaxis()->SetTitle("counts");
394   
395   fHistTOFvsEcells = new TH2F("fHistTOFvsEcells","TOF vs Energy of cells", TOFbins, TOFlow, TOFup, Ebins, Elow, Eup);
396   fHistTOFvsEcells->GetXaxis()->SetTitle("TOF [s]");
397   fHistTOFvsEcells->GetYaxis()->SetTitle("Energy [GeV]");
398         
399   fHistTOFvsEcellsC = new TH2F("fHistTOFvsEcellsC","TOF vs Energy of cells (corrected)", TOFbins, TOFlow, TOFup, Ebins, Elow, Eup);
400   fHistTOFvsEcellsC->GetXaxis()->SetTitle("TOF");
401   fHistTOFvsEcellsC->GetYaxis()->SetTitle("E [GeV]");
402   
403   fHistEmaxCellvsAmpFastOR = new TH2F("fHistEmaxCellvsAmpFastOR","Maximum Cell Energy vs Maximum FastOR Amplitude", L0Ampbins, L0Amplow, L0Ampup, Ebins, Elow, Eup);
404   fHistEmaxCellvsAmpFastOR->GetXaxis()->SetTitle("Maximum FastOR Amplitude");
405   fHistEmaxCellvsAmpFastOR->GetYaxis()->SetTitle("Maximum Cell Energy [GeV]");
406   
407   fOutput->Add(fHistEclus);
408   fOutput->Add(fHistEmaxClus);
409   fOutput->Add(fHistEtavsPhiMaxClus);
410   fOutput->Add(fHistEtavsEmaxClus);
411   fOutput->Add(fHistPhivsEmaxClus);
412   fOutput->Add(fHistTOFvsEclus);
413   fOutput->Add(fHistTOFvsEclusC);
414   fOutput->Add(fHistNcellsvsEclus);
415   fOutput->Add(fHistAmpTClus);
416   fOutput->Add(fHistAmpMaxTClus);
417   fOutput->Add(fHistEtavsPhiMaxTClus);
418   fOutput->Add(fHistEmaxClusvsAmpMaxTClus);
419   fOutput->Add(fHistEmaxClusvsAmpMatchedTClus);
420   fOutput->Add(fHistEmaxClusNotMatchingTClus);
421   fOutput->Add(fHistEtavsPhiMaxClusNotMatchingTClus);
422   fOutput->Add(fHistEmatchedClusvsAmpMaxTClus);
423   fOutput->Add(fHistAmpMaxTClusNotMatchingClus);
424   fOutput->Add(fHistEtavsPhiMaxTClusNotMatchingClus);
425   fOutput->Add(fHistIdxMaxClusvsIdxMaxTClus);
426   fOutput->Add(fHistPhiMaxClusvsPhiMaxTClus);
427   fOutput->Add(fHistEtaMaxClusvsEtaMaxTClus);
428   fOutput->Add(fHistTOFmaxClusvsTimeMaxTClus);
429   fOutput->Add(fHistEmatchedClusvsAmpMatchedTClus);
430   fOutput->Add(fHistEmatchedClus);
431   fOutput->Add(fHistEmaxMatchedClus);
432   fOutput->Add(fHistAmpL1TimeSum);
433   fOutput->Add(fHistAmpMaxL1TimeSum);
434   fOutput->Add(fHistAmpMaxL1TimeSumVScent);
435   fOutput->Add(fHistAmpFastORvsAmpL1TimeSum);
436   fOutput->Add(fHistAmpFastOR);
437   fOutput->Add(fHistAmpMaxFastOR);
438   fOutput->Add(fHistTimeFastOR);
439   fOutput->Add(fHistEtavsPhiFastOR);
440   fOutput->Add(fHistEtavsPhiMaxFastOR);
441   fOutput->Add(fHistTimeDispFastOR);
442   fOutput->Add(fHistTimevsL0TimeFastOR);
443   fOutput->Add(fHistNtimesFastOR);
444   fOutput->Add(fHistEcells);
445   fOutput->Add(fHistEmaxCell);
446   fOutput->Add(fHistTOFvsEcells);
447   fOutput->Add(fHistTOFvsEcellsC);
448   fOutput->Add(fHistEmaxCellvsAmpFastOR);
449         
450   PostData(1, fOutput); 
451 }
452
453 //________________________________________________________________________
454 void AliAnalysisTaskSATR::Init() 
455 {
456   fGeom = AliEMCALGeometry::GetInstance("EMCAL_COMPLETEV1");
457   if (!fGeom) {
458     AliError("Couldn't get geometry. Returning...");
459     return;
460   }
461
462   if (fRun <= 0)
463     return;
464   
465   AliCDBManager *cdb = AliCDBManager::Instance();
466   
467   if (!fPedestal && fLoadPed) {
468     if (!cdb->IsDefaultStorageSet() && !fOCDBpath.IsNull())
469       cdb->SetDefaultStorage(fOCDBpath);
470     cdb->SetRun(fRun);
471     AliCDBEntry *entry = static_cast<AliCDBEntry*>(AliCDBManager::Instance()->Get("EMCAL/Calib/Pedestals"));
472     if (entry) 
473       fPedestal =  static_cast<AliCaloCalibPedestal*>(entry->GetObject());
474   }
475 }
476
477 //________________________________________________________________________
478 void AliAnalysisTaskSATR::UserExec(Option_t *) 
479 {
480   // Main loop
481   // Called for each event
482   
483   // Create pointer to reconstructed event
484   AliVEvent *event = InputEvent();
485   if (!event){ 
486     AliError("ERROR: Could not retrieve event"); 
487     return; 
488   }
489   
490   fRun = event->GetRunNumber();
491
492   Float_t cent = 99;
493
494   AliCentrality *centobj = InputEvent()->GetCentrality();
495   if (centobj) {
496     cent = centobj->GetCentralityPercentile("V0M");
497   }
498   
499   Init();
500   
501   TClonesArray *caloClusters = dynamic_cast<TClonesArray*>(event->FindListObject(fCaloClustersName));  
502
503   if (!caloClusters){
504     AliError("Troubles trying to load clusters!");
505     caloClusters = dynamic_cast<TClonesArray*>(event->FindListObject("caloClusters"));
506   }
507  
508   if(!caloClusters){
509     AliError("Cannot get the CaloClusters!");
510     return;
511   }
512   
513   Int_t NumberOfCaloClusters = caloClusters->GetEntriesFast();
514   
515   TClonesArray *triggerClusters = dynamic_cast<TClonesArray*>(event->FindListObject(fTriggerClustersName));  
516  
517   Int_t NumberOfTriggerClusters = 0;
518   
519   if(triggerClusters){
520     NumberOfTriggerClusters = triggerClusters->GetEntriesFast();
521   }
522   
523   Float_t clusterTime_corr;
524   if(event->GetBunchCrossNumber() % 4 < 2 )
525     clusterTime_corr = -0.0000001;
526   else 
527     clusterTime_corr = 0;
528   
529   Float_t maxL0amp = -1;
530   
531   for (Int_t i = 0; i < NumberOfTriggerClusters; i++) {
532     AliVCluster* trgCluster = static_cast<AliVCluster*>(triggerClusters->At(i));
533     if (!trgCluster)
534       continue;
535     Float_t L0amp = trgCluster->E();
536
537     if (fCheckDeadClusters && trgCluster->GetDistanceToBadChannel() < 3)
538       continue;
539     
540     if (fTimeCutOn && (trgCluster->GetTOF() < fMinL0Time || trgCluster->GetTOF() > fMaxL0Time))
541       continue;
542     
543     if (maxL0amp < L0amp)
544       maxL0amp = L0amp;
545     } // trigger cluster loop
546   
547   if (maxL0amp > fMaxCutL0Amp || maxL0amp < fMinCutL0Amp) {
548     AliWarning(Form("Event skipped because maximum trigger cluster amplitude = %f out of limits (%f, %f)", maxL0amp, fMinCutL0Amp, fMaxCutL0Amp));
549     return; 
550   }
551   
552   Float_t maxe = -1;
553   
554   for (Int_t iClusters = 0; iClusters < NumberOfCaloClusters; iClusters++) {
555     AliVCluster* cluster = dynamic_cast<AliVCluster*>(caloClusters->At(iClusters));
556     if (!cluster){
557       AliError(Form("Could not receive cluster %d", iClusters));
558       continue;
559     }  
560     
561     if (!(cluster->IsEMCAL())) 
562       continue;
563     
564     if (fCheckDeadClusters && cluster->GetDistanceToBadChannel() < 3) {
565       continue;
566     }
567     
568     if (fTimeCutOn && (cluster->GetTOF() + clusterTime_corr < fMinClusTime || cluster->GetTOF() + clusterTime_corr > fMaxClusTime))
569       continue;
570     
571     if (maxe < cluster->E()){
572       maxe = cluster->E();
573     }
574   } //cluster loop 
575     
576   if (maxe > fMaxCutClusEnergy || maxL0amp < fMinCutClusEnergy) {
577     AliWarning(Form("Event skipped because maximum cluster energy = %f out of limits (%f, %f)", maxe, fMaxCutClusEnergy, fMinCutClusEnergy));
578     return;
579   }
580   
581   AliVCaloTrigger *triggers = event->GetCaloTrigger("EMCAL");
582   
583   Float_t maxL0FastORamp = -1;
584   Float_t maxL0FastORtime = -1;
585   
586   if (triggers && triggers->GetEntries() > 0) {
587     
588     triggers->Reset();
589     Float_t L0FastORamp = 0;
590     
591     while (triggers->Next()) {
592       
593       triggers->GetAmplitude(L0FastORamp);
594       
595       if (L0FastORamp < 0)
596         continue;
597       
598       if (maxL0FastORamp < L0FastORamp) {
599         
600         Int_t ntimes = 0;
601         
602         triggers->GetNL0Times(ntimes);
603         
604         if (fTimeCutOn && ntimes < 1) {
605           AliWarning(Form("FastOR removed from analysis because did not contribute to L0 trigger (no L0 time information)"));
606           continue;
607         }
608         
609         maxL0FastORamp = L0FastORamp;
610         
611         Int_t trgtimes[25];
612         triggers->GetL0Times(trgtimes);
613         
614         maxL0FastORtime = trgtimes[0];
615         
616       }
617     }
618   }
619   
620   maxL0FastORamp = -1;
621   Int_t maxL0FastORrow = -1;
622   Int_t maxL0FastORcol = -1;
623   Float_t maxL1amp = -1;
624   Float_t TimeBinAmp = 0;
625   
626   if (triggers && triggers->GetEntries() > 0) {
627     
628     triggers->Reset();
629     Float_t L0FastORamp = 0;
630     Int_t L1amp = 0;
631     
632     while (triggers->Next()) {
633       
634       triggers->GetAmplitude(L0FastORamp);
635       
636       if (L0FastORamp < 0)
637         continue;
638       
639       Int_t ntimes = 0;
640       
641       triggers->GetNL0Times(ntimes);
642       
643       triggers->GetTime(TimeBinAmp);
644       
645       fHistNtimesFastOR->Fill(ntimes);
646       
647       if (ntimes > 0) {
648         Int_t trgtimes[25];
649         triggers->GetL0Times(trgtimes);
650         
651         Int_t mintime = trgtimes[0];
652         Int_t maxtime = trgtimes[0];
653         
654         for (Int_t i = 0; i < ntimes; ++i) {
655           if (trgtimes[i] < mintime)
656             mintime = trgtimes[i];
657           if (maxtime < trgtimes[i])
658             maxtime = trgtimes[i];
659           
660           fHistTimeFastOR->Fill(trgtimes[i]);
661           fHistTimevsL0TimeFastOR->Fill(TimeBinAmp, trgtimes[i]);
662           fHistTimeDispFastOR->Fill(trgtimes[i] - maxL0FastORtime);
663         }
664         
665         if (fTimeCutOn && ((fMinL0Time > mintime) || (fMaxL0Time < maxtime)))
666           continue;
667       }
668       
669       if (fTimeCutOn && ntimes < 1)
670         continue;
671       
672       Int_t gCol=0, gRow=0;
673       triggers->GetPosition(gCol, gRow);
674       
675       Int_t find = -1;
676       fGeom->GetAbsFastORIndexFromPositionInEMCAL(gCol,gRow,find);
677       
678       if (find < 0)
679         continue;
680       
681       Int_t cidx[4] = {-1};
682       Bool_t ret = fGeom->GetCellIndexFromFastORIndex(find, cidx);
683       
684       if (!ret)
685         continue;
686
687       Int_t nSupMod = 0, nModule = 0, nIphi = 0, nIeta = 0, iphi = 0, ieta = 0;      
688
689       Bool_t deadCluster = kFALSE;
690       
691       if (fCheckDeadClusters && fPedestal) {
692         for (Int_t i = 0; i < 4; i++){
693           fGeom->GetCellIndex(cidx[i], nSupMod, nModule, nIphi, nIeta);
694           fGeom->GetCellPhiEtaIndexInSModule(nSupMod, nModule, nIphi, nIeta, iphi, ieta);
695           
696           Double_t d = fPedestal->GetDeadMap(nSupMod)->GetBinContent(ieta,iphi);
697           if (d == AliCaloCalibPedestal::kDead || d == AliCaloCalibPedestal::kHot){
698             AliWarning(Form("Dead/hot FastOR removed from analysis"));
699             deadCluster = kTRUE;
700           }
701         }
702       }
703       
704       if (deadCluster)
705         continue;
706
707       fHistAmpFastOR->Fill(L0FastORamp);
708       fHistEtavsPhiFastOR->Fill(gCol, gRow);
709       
710       if (maxL0FastORamp < L0FastORamp) {
711         maxL0FastORamp = L0FastORamp;
712         maxL0FastORcol = gCol;
713         maxL0FastORrow = gRow;
714       }
715       
716       triggers->GetL1TimeSum(L1amp);
717       fHistAmpL1TimeSum->Fill(L1amp);
718       fHistAmpFastORvsAmpL1TimeSum->Fill(L0FastORamp,L1amp);
719       
720       if (maxL1amp < L1amp) 
721         maxL1amp = L1amp;
722   
723     }
724   }
725   
726   if (maxL0FastORamp > -1) {
727     fHistAmpMaxFastOR->Fill(maxL0FastORamp); 
728     fHistEtavsPhiMaxFastOR->Fill(maxL0FastORcol, maxL0FastORrow);
729   }
730   
731   if (maxL1amp > -1) {
732     fHistAmpMaxL1TimeSum->Fill(maxL1amp);
733     fHistAmpMaxL1TimeSumVScent->Fill(cent, maxL1amp);
734   }
735   
736   maxL0amp = 0;
737   Float_t EmatchMaxL0 = -1;
738   Int_t maxL0index = -1;
739   Float_t EmatchL0 = -1;
740   AliVCluster* maxTcluster = 0;
741  
742   for (Int_t i = 0; i < NumberOfTriggerClusters; i++){
743     AliVCluster* trgCluster = dynamic_cast<AliVCluster*>(triggerClusters->At(i));
744     if (!trgCluster)
745       continue;
746     if (fCheckDeadClusters && trgCluster->GetDistanceToBadChannel() < 3) {
747       AliWarning(Form("Trigger cluster %d removed from analysis because distance to bad channel = %f < 3", i, trgCluster->GetDistanceToBadChannel()));
748       continue;
749     }
750     
751     if (fTimeCutOn && (trgCluster->GetTOF() < fMinL0Time || trgCluster->GetTOF() > fMaxL0Time)) {
752       AliWarning(Form("Trigger cluster %d removed from analysis because out of time %f, limits (%d, %d)", i, trgCluster->GetTOF(), fMinL0Time, fMaxL0Time));
753       continue;
754     }
755     Float_t L0amp = trgCluster->E();
756     
757     fHistAmpTClus->Fill(L0amp);
758       
759     AliVCluster* matchedcluster = GetClusterFromId(caloClusters, trgCluster->GetID());
760     
761     if (matchedcluster) {
762       EmatchL0 = matchedcluster->E();
763
764       fHistEmatchedClus->Fill(EmatchL0);
765       fHistEmatchedClusvsAmpMatchedTClus->Fill(L0amp, EmatchL0);
766     }
767     
768     if (maxL0amp < L0amp) {
769       maxL0amp = L0amp;
770       maxL0index = trgCluster->GetID();
771       EmatchMaxL0 = EmatchL0;
772       maxTcluster = trgCluster;
773     }
774   }
775   
776   Int_t etaMaxTClus = 0;
777   Int_t phiMaxTClus = 0;
778   
779   if (maxTcluster != 0) {
780     fHistAmpMaxTClus->Fill(maxL0amp);
781     fHistEmaxMatchedClus->Fill(EmatchMaxL0);
782     
783     if (fTriggerClusterizer) {
784       Int_t nEtaDigitsSupMod = fGeom->GetNEta() * fGeom->GetNETAdiv(); // always 48?;
785       Int_t nPhiDigitsSupMod = fGeom->GetNPhi() * fGeom->GetNPHIdiv(); // always 24?;
786       
787       Int_t nTRUPhi = 1;
788       Int_t nTRUEta = 1;
789       
790       Int_t nEtaDigits = nEtaDigitsSupMod * fGeom->GetNumberOfSuperModules() / fGeom->GetNPhiSuperModule();
791       Int_t nPhiDigits = nPhiDigitsSupMod * fGeom->GetNPhiSuperModule();    
792       
793       if (fTriggerClusterizer->GetTRUShift()) {
794         nTRUPhi = fGeom->GetNPhiSuperModule() * 3;
795         nTRUEta = fGeom->GetNumberOfSuperModules() / fGeom->GetNPhiSuperModule();
796         nEtaDigits /= nTRUEta;
797         nPhiDigits /= nTRUPhi;
798       }
799       
800       Int_t nClusEtaNoShift = nEtaDigits / fTriggerClusterizer->GetnEta();
801       Int_t nClusPhiNoShift = nPhiDigits / fTriggerClusterizer->GetnPhi();
802       
803       Int_t nClusters =  nClusEtaNoShift * nClusPhiNoShift * nTRUEta * nTRUPhi;
804        
805       etaMaxTClus = (maxTcluster->GetID() % nClusters) / (nClusPhiNoShift * nTRUPhi);
806       phiMaxTClus = (maxTcluster->GetID() % nClusters) % (nClusPhiNoShift * nTRUPhi);
807       
808       fHistEtavsPhiMaxTClus->Fill(etaMaxTClus, phiMaxTClus);
809     }
810     
811     AliVCluster *matchedClus = GetClusterFromId(caloClusters, maxTcluster->GetID());
812     
813     if (!matchedClus) {
814       fHistAmpMaxTClusNotMatchingClus->Fill(maxL0amp);
815       if (fClusterizer)
816         fHistEtavsPhiMaxTClusNotMatchingClus->Fill(etaMaxTClus, phiMaxTClus);
817     }
818     else {
819       fHistEmatchedClusvsAmpMaxTClus->Fill(maxL0amp, matchedClus->E());
820     }
821   }
822     
823   maxe = 0;
824   Float_t maxClusterTime = 0;
825         AliVCluster* maxcluster = 0;
826   Int_t maxClusId = -1;
827   Int_t maxiCluster = -1;
828   
829   for (Int_t iClusters = 0; iClusters < NumberOfCaloClusters; iClusters++) {
830     AliVCluster* cluster = dynamic_cast<AliVCluster*>(caloClusters->At(iClusters));
831     if (!cluster){
832       AliError(Form("Could not receive cluster %d", iClusters));
833       continue;
834     }  
835     
836     if (!(cluster->IsEMCAL())) 
837       continue;
838     
839     if (fCheckDeadClusters && cluster->GetDistanceToBadChannel() < 3) {
840       AliWarning(Form("Cluster %d removed from analysis because distance to bad channel = %f < 3", iClusters, cluster->GetDistanceToBadChannel()));
841       continue;
842     }
843     
844     if (fTimeCutOn && (cluster->GetTOF() < fMinL0Time || cluster->GetTOF() > fMaxL0Time)) {
845       AliWarning(Form("Cluster %d removed from analysis because out of time %f, limits (%d, %d)", iClusters, cluster->GetTOF(), fMinL0Time, fMaxL0Time));
846       continue;
847     }
848     
849     if (maxe < cluster->E()){
850       maxe = cluster->E();
851       maxcluster = cluster;
852       maxClusterTime = cluster->GetTOF() + clusterTime_corr;
853       maxClusId = cluster->GetID();
854       maxiCluster = iClusters;
855     }
856     
857     fHistEclus->Fill(cluster->E());
858     
859     fHistTOFvsEclus->Fill(cluster->GetTOF(),cluster->E());
860     
861     fHistTOFvsEclusC->Fill(cluster->GetTOF() + clusterTime_corr,cluster->E());
862     
863     fHistNcellsvsEclus->Fill(cluster->GetNCells(),cluster->E());
864     
865   } //cluster loop 
866   
867   Int_t etaMaxClus = 0;
868   Int_t phiMaxClus = 0;
869   
870   if (maxcluster != 0) {
871     
872     fHistEmaxClus->Fill(maxe);
873     
874     if (fClusterizer) {
875       Int_t nEtaDigitsSupMod = fGeom->GetNEta() * fGeom->GetNETAdiv(); // always 48?;
876       Int_t nPhiDigitsSupMod = fGeom->GetNPhi() * fGeom->GetNPHIdiv(); // always 24?;
877       
878       Int_t nTRUPhi = 1;
879       Int_t nTRUEta = 1;
880       
881       Int_t nEtaDigits = nEtaDigitsSupMod * fGeom->GetNumberOfSuperModules() / fGeom->GetNPhiSuperModule();
882       Int_t nPhiDigits = nPhiDigitsSupMod * fGeom->GetNPhiSuperModule();    
883       
884       if (fClusterizer->GetTRUShift()) {
885         nTRUPhi = fGeom->GetNPhiSuperModule() * 3;
886         nTRUEta = fGeom->GetNumberOfSuperModules() / fGeom->GetNPhiSuperModule();
887         nEtaDigits /= nTRUEta;
888         nPhiDigits /= nTRUPhi;
889       }
890       
891       Int_t nClusEtaNoShift = nEtaDigits / fClusterizer->GetnEta();
892       Int_t nClusPhiNoShift = nPhiDigits / fClusterizer->GetnPhi();
893       
894       Int_t nClusters =  nClusEtaNoShift * nClusPhiNoShift * nTRUEta * nTRUPhi;
895       
896       etaMaxClus = (maxcluster->GetID() % nClusters) / (nClusPhiNoShift * nTRUPhi);
897       phiMaxClus = (maxcluster->GetID() % nClusters) % (nClusPhiNoShift * nTRUPhi);
898       
899       fHistEtavsPhiMaxClus->Fill(etaMaxClus, phiMaxClus);
900       fHistEtavsEmaxClus->Fill(etaMaxClus, maxcluster->E());
901       fHistPhivsEmaxClus->Fill(phiMaxClus, maxcluster->E());
902     }
903     
904     AliVCluster *matchedTClus = GetClusterFromId(triggerClusters, maxcluster->GetID());
905     
906     if (!matchedTClus || (fTimeCutOn && (matchedTClus->GetTOF() < fMinL0Time || matchedTClus->GetTOF() > fMaxL0Time))) {
907       fHistEmaxClusNotMatchingTClus->Fill(maxe);
908       if (fClusterizer)
909         fHistEtavsPhiMaxClusNotMatchingTClus->Fill(etaMaxClus, phiMaxClus);
910     }
911     else {
912       fHistEmaxClusvsAmpMatchedTClus->Fill(matchedTClus->E(), maxe);
913     }
914   }
915         
916   if (maxTcluster != 0 && maxcluster != 0) {    
917     fHistEmaxClusvsAmpMaxTClus->Fill(maxL0amp, maxe);    
918     fHistTOFmaxClusvsTimeMaxTClus->Fill(maxTcluster->GetTOF(), maxcluster->GetTOF() + clusterTime_corr);
919     fHistIdxMaxClusvsIdxMaxTClus->Fill(maxL0index, maxClusId);
920     fHistPhiMaxClusvsPhiMaxTClus->Fill(phiMaxTClus, phiMaxClus);
921     fHistEtaMaxClusvsEtaMaxTClus->Fill(etaMaxTClus, etaMaxClus);
922   }
923   
924   AliVCaloCells *cells = event->GetEMCALCells();
925   Int_t icells;
926   Float_t maxCellEnergy = -1;
927   for (icells = 0; icells < cells->GetNumberOfCells(); icells++) {
928     if (fTimeCutOn && (cells->GetTime(icells) + clusterTime_corr < fMinClusTime || cells->GetTime(icells) + clusterTime_corr > fMaxClusTime))
929       continue;
930     
931     fHistEcells->Fill(cells->GetAmplitude(icells));
932     
933     fHistTOFvsEcells->Fill(cells->GetTime(icells), cells->GetAmplitude(icells));
934     fHistTOFvsEcellsC->Fill(cells->GetTime(icells) + clusterTime_corr, cells->GetAmplitude(icells));
935     //cout << cells->GetTime(icells) << endl;
936     
937     if (maxCellEnergy < cells->GetAmplitude(icells)) 
938       maxCellEnergy = cells->GetAmplitude(icells);
939     
940   }
941   
942   if (maxCellEnergy > -1) {
943     fHistEmaxCell->Fill(maxCellEnergy);
944   }
945   
946   if (maxL0FastORamp > -1 && maxCellEnergy > -1) {
947     fHistEmaxCellvsAmpFastOR->Fill(maxL0FastORamp, maxCellEnergy);
948   }
949   
950   // information for this iteration of the UserExec in the container
951   PostData(1, fOutput);
952   
953 }
954
955
956 //________________________________________________________________________
957 void AliAnalysisTaskSATR::Terminate(Option_t *) 
958 {
959
960 }
961
962 AliVCluster* AliAnalysisTaskSATR::GetClusterFromId(TClonesArray *caloClusters, Int_t id)
963 {
964   for (Int_t iClusters = 0; iClusters < caloClusters->GetEntriesFast(); iClusters++){
965     AliVCluster* cluster = dynamic_cast<AliVCluster*>(caloClusters->At(iClusters));
966     if (!cluster){
967       AliError(Form("Could not receive cluster %d", iClusters));
968       continue;
969     }  
970     
971     if (!(cluster->IsEMCAL())) 
972       continue;
973     
974     if (id == cluster->GetID())
975       return cluster;
976   } //cluster loop 
977   
978   return 0;
979 }