]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/AliZDCQADataMakerRec.cxx
Fiexs needed for Pb-Pb reconstruction
[u/mrichter/AliRoot.git] / ZDC / AliZDCQADataMakerRec.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 ///////////////////////////////////////////////////////////////////////
17 //                                                                   //
18 //  Produces the data needed to calculate the ZDC quality assurance. //
19 //  QA objects are 1 & 2 Dimensional histograms.                     //
20 //  author: C. Oppedisano                                            //
21 //                                                                   //
22 ///////////////////////////////////////////////////////////////////////
23
24 // --- ROOT system ---
25 #include <TClonesArray.h>
26 #include <TFile.h> 
27 #include <TH1F.h> 
28 #include <TH2F.h>
29 #include <TLine.h>
30 #include <TProfile.h>
31 #include <Riostream.h>
32 // --- Standard library ---
33
34 // --- AliRoot header files ---
35 #include "AliLog.h"
36 #include "AliQAChecker.h"
37 #include "AliZDCReco.h"
38 #include "AliRawReader.h"
39 #include "AliZDCQADataMakerRec.h"
40 #include "AliZDCPedestals.h"
41 #include "AliZDCRawStream.h"
42 #include "AliZDCDigit.h"
43 #include "AliESDZDC.h"
44 #include "AliESDEvent.h"
45
46 ClassImp(AliZDCQADataMakerRec)
47            
48 //____________________________________________________________________________ 
49   AliZDCQADataMakerRec::AliZDCQADataMakerRec() : 
50   AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kZDC), "ZDC Quality Assurance Data Maker"), 
51   fPedCalibData(0x0)
52 {
53   // ctor
54 }
55
56 //____________________________________________________________________________ 
57 AliZDCQADataMakerRec::AliZDCQADataMakerRec(const AliZDCQADataMakerRec& qadm) :
58   AliQADataMakerRec(),      
59   fPedCalibData(qadm.fPedCalibData)
60
61 {
62   //copy ctor 
63   SetName((const char*)qadm.GetName()); 
64   SetTitle((const char*)qadm.GetTitle()); 
65 }
66
67 //__________________________________________________________________
68 AliZDCQADataMakerRec& AliZDCQADataMakerRec::operator = (const AliZDCQADataMakerRec& qadm )
69 {
70   // Equal operator.
71   this->~AliZDCQADataMakerRec();
72   new(this) AliZDCQADataMakerRec(qadm);
73   return *this;
74 }
75
76 //____________________________________________________________________________ 
77 AliZDCQADataMakerRec::~AliZDCQADataMakerRec()
78 {
79   if(fPedCalibData && !(AliCDBManager::Instance()->GetCacheFlag())){
80     delete fPedCalibData;
81     fPedCalibData=0;
82   } 
83 }
84
85 //____________________________________________________________________________ 
86 AliZDCPedestals* AliZDCQADataMakerRec::GetPedCalibData() const
87 {
88
89   // Retrieving pedestal calibration object from OCDB
90   AliCDBEntry  *entry = AliCDBManager::Instance()->Get("ZDC/Calib/Pedestals");
91   if(!entry) AliWarning("No calibration data loaded!");  
92
93   AliZDCPedestals *calibdata = (AliZDCPedestals*)  (entry->GetObject());
94   if(!calibdata) AliFatal("Wrong calibration object in calibration  file!");
95
96   return calibdata;
97
98 }
99
100 //____________________________________________________________________________ 
101 void AliZDCQADataMakerRec::InitDigits()
102 {
103   // create Digits histograms in Digits subdir
104   //
105   const Bool_t expert   = kTRUE ; 
106   const Bool_t image    = kTRUE ; 
107   
108   // ------------------- HIGH GAIN CHAIN ---------------------------
109   TH1F * hDigZNCTot = new TH1F("hDigZNCTot", "Signal in ZNC;Amplitude [ADC counts];Counts", 100, 0., 6000.);
110   TH1F * hDigZNATot = new TH1F("hDigZNATot", "Signal in ZNA;Amplitude [ADC counts];Counts", 100, 0., 6000.);
111   TH1F * hDigZPCTot = new TH1F("hDigZPCTot", "Signal in ZPC;Amplitude [ADC counts];Counts", 100, 0., 6000.);
112   TH1F * hDigZPATot = new TH1F("hDigZPATot", "Signal in ZPA;Amplitude [ADC counts];Counts", 100, 0., 6000.);
113   Add2DigitsList(hDigZNCTot, 0, !expert, image);
114   Add2DigitsList(hDigZNATot, 1, !expert, image);
115   Add2DigitsList(hDigZPCTot, 2, !expert, image);
116   Add2DigitsList(hDigZPATot, 3, !expert, image);
117   //
118   TH1F * hDigSumQZNC = new TH1F("hDigSumQZNC", "Signal in 4 ZNC PMQ;Amplitude [ADC counts];Counts",100, 0., 4000.);
119   TH1F * hDigSumQZNA = new TH1F("hDigSumQZNA", "Signal in 4 ZNA PMQ;Amplitude [ADC counts];Counts",100, 0., 4000.);
120   TH1F * hDigSumQZPC = new TH1F("hDigSumQZPC", "Signal in 4 ZPC PMQ;Amplitude [ADC counts];Counts",100, 0., 4000.);
121   TH1F * hDigSumQZPA = new TH1F("hDigSumQZPA", "Signal in 4 ZPA PMQ;Amplitude [ADC counts];Counts",100, 0., 4000.);
122   Add2DigitsList(hDigSumQZNC, 4, expert, !image);
123   Add2DigitsList(hDigSumQZNA, 5, expert, !image);
124   Add2DigitsList(hDigSumQZPC, 6, expert, !image);
125   Add2DigitsList(hDigSumQZPA, 7, expert, !image);
126   //
127   TH1F * hDigPMCZNC = new TH1F("hDigPMCZNC", "Signal in ZNC PMC;Amplitude [ADC counts];Counts",100, 0., 4000.);
128   TH1F * hDigPMCZNA = new TH1F("hDigPMCZNA", "Signal in ZNA PMC;Amplitude [ADC counts];Counts",100, 0., 4000.);
129   TH1F * hDigPMCZPC = new TH1F("hDigPMCZPC", "Signal in ZPC PMC;Amplitude [ADC counts];Counts",100, 0., 4000.);
130   TH1F * hDigPMCZPA = new TH1F("hDigPMCZPA", "Signal in ZPA PMC;Amplitude [ADC counts];Counts",100, 0., 4000.);
131   Add2DigitsList(hDigPMCZNC, 8, expert, !image);
132   Add2DigitsList(hDigPMCZNA, 9, expert, !image);
133   Add2DigitsList(hDigPMCZPC, 10, expert, !image);
134   Add2DigitsList(hDigPMCZPA, 11, expert, !image);
135   // 
136   // ------------------- LOW GAIN CHAIN ---------------------------
137   TH1F * hDigZNCTotlg = new TH1F("hDigZNCTotlg", "Digit lg signal in ZNC", 100, 0., 6000.);
138   TH1F * hDigZNATotlg = new TH1F("hDigZNATotlg", "Digit lg signal in ZNA", 100, 0., 6000.);
139   TH1F * hDigZPCTotlg = new TH1F("hDigZPCTotlg", "Digit lg signal in ZPC", 100, 0., 6000.);
140   TH1F * hDigZPATotlg = new TH1F("hDigZPATotlg", "Digit lg signal in ZPA", 100, 0., 6000.);
141   Add2DigitsList(hDigZNCTotlg, 12, expert, !image);
142   Add2DigitsList(hDigZNATotlg, 13, expert, !image);
143   Add2DigitsList(hDigZPCTotlg, 14, expert, !image);
144   Add2DigitsList(hDigZPATotlg, 15, expert, !image);
145   //
146   TH1F * hDigSumQZNClg = new TH1F("hDigSumQZNClg", "Signal in 4 ZNC PMQlg",100, 0., 4000.);
147   TH1F * hDigSumQZNAlg = new TH1F("hDigSumQZNAlg", "Signal in 4 ZNA PMQlg",100, 0., 4000.);
148   TH1F * hDigSumQZPClg = new TH1F("hDigSumQZPClg", "Signal in 4 ZPC PMQlg",100, 0., 4000.);
149   TH1F * hDigSumQZPAlg = new TH1F("hDigSumQZPAlg", "Signal in 4 ZPA PMQlg",100, 0., 4000.);
150   Add2DigitsList(hDigSumQZNClg, 16, expert, !image);
151   Add2DigitsList(hDigSumQZNAlg, 17, expert, !image);
152   Add2DigitsList(hDigSumQZPClg, 18, expert, !image);
153   Add2DigitsList(hDigSumQZPAlg, 19, expert, !image);
154   //
155   TH1F * hDigPMCZNClg = new TH1F("hDigPMCZNClg", "Signal in ZNC PMClg",100, 0., 4000.);
156   TH1F * hDigPMCZNAlg = new TH1F("hDigPMCZNAlg", "Signal in ZNA PMClg",100, 0., 4000.);
157   TH1F * hDigPMCZPClg = new TH1F("hDigPMCZPClg", "Signal in ZPC PMClg",100, 0., 4000.);
158   TH1F * hDigPMCZPAlg = new TH1F("hDigPMCZPAlg", "Signal in ZPA PMClg",100, 0., 4000.);
159   Add2DigitsList(hDigPMCZNClg, 20, expert, !image);
160   Add2DigitsList(hDigPMCZNAlg, 21, expert, !image);
161   Add2DigitsList(hDigPMCZPClg, 22, expert, !image);
162   Add2DigitsList(hDigPMCZPAlg, 23, expert, !image);
163
164 }
165
166 //____________________________________________________________________________
167 void AliZDCQADataMakerRec::InitRaws()
168 {
169   // create Digits histograms in Digits subdir
170   const Bool_t expert   = kTRUE ; 
171   const Bool_t image    = kTRUE ; 
172
173   TH1F * hZNCSpectrum = new TH1F("hZNCSpectrum","ZNC spectrum;Amplitude [ADC counts];Counts",100,0.,1200.);
174   TH1F * hZPCSpectrum = new TH1F("hZPCSpectrum","ZPC spectrum;Amplitude [ADC counts];Counts",100,0.,1200.);
175   TH1F * hZNASpectrum = new TH1F("hZNASpectrum","ZNA spectrum;Amplitude [ADC counts];Counts",100,0.,1200.);
176   TH1F * hZPASpectrum = new TH1F("hZPASpectrum","ZPA spectrum;Amplitude [ADC counts];Counts",100,0.,1200.);
177   TH1F * hZEM1Spectrum = new TH1F("hZEM1Spectrum","ZEM1 spectrum;Amplitude [ADC counts];Counts",100,8., 1208.);
178   TH1F * hZEM2Spectrum = new TH1F("hZEM2Spectrum","ZEM2 spectrum;Amplitude [ADC counts];Counts",100,8., 1208.);
179   Add2RawsList(hZNCSpectrum, 0, expert, !image);
180   Add2RawsList(hZNASpectrum, 1, expert, !image);
181   Add2RawsList(hZPCSpectrum, 2, expert, !image);
182   Add2RawsList(hZPASpectrum, 3, expert, !image);
183   Add2RawsList(hZEM1Spectrum, 4, expert, !image);
184   Add2RawsList(hZEM2Spectrum, 5, expert, !image);
185   //
186   TH2F * hZNCpmCvsPMq = new TH2F("hZNCpmCvsPMq", "ZNC;PMC [ADC counts];Sum(PMQ) [ADC counts]",50,8.,1208.,50,8., 1208.);
187   TH2F * hZPCpmCvsPMq = new TH2F("hZPCpmCvsPMq", "ZPC;PMC [ADC counts];Sum(PMQ) [ADC counts]",50,8.,1208.,50,8., 1208.);
188   TH2F * hZNApmCvsPMq = new TH2F("hZNApmCvsPMq", "ZNA;PMC [ADC counts];Sum(PMQ) [ADC counts]",50,8.,1208.,50,8., 1208.);
189   TH2F * hZPApmCvsPMq = new TH2F("hZPApmCvsPMq", "ZPA;PMC [ADC counts];Sum(PMQ) [ADC counts]",50,8.,1208.,50,8., 1208.);
190   Add2RawsList(hZNCpmCvsPMq, 6, expert, !image);
191   Add2RawsList(hZNApmCvsPMq, 7, expert, !image);
192   Add2RawsList(hZPCpmCvsPMq, 8, expert, !image);
193   Add2RawsList(hZPApmCvsPMq, 9, expert, !image);
194     
195   TH1F * hRawPMCZNC = new TH1F("hRawPMCZNC", "Raw ZNC PMC;Amplitude [ADC counts];Counts",100, 8., 1208.);
196   TH1F * hRawPMCZNA = new TH1F("hRawPMCZNA", "Raw ZNA PMC;Amplitude [ADC counts];Counts",100, 8., 1208.);
197   TH1F * hRawPMCZPC = new TH1F("hRawPMCZPC", "Raw ZPC PMC;Amplitude [ADC counts];Counts",100, 8., 1208.);
198   TH1F * hRawPMCZPA = new TH1F("hRawPMCZPA", "Raw ZPA PMC;Amplitude [ADC counts];Counts",100, 8., 1208.);
199   Add2RawsList(hRawPMCZNC, 10, expert, !image);
200   Add2RawsList(hRawPMCZNA, 11, expert, !image);
201   Add2RawsList(hRawPMCZPC, 12, expert, !image);
202   Add2RawsList(hRawPMCZPA, 13, expert, !image);
203   TH1F * hRawSumQZNC = new TH1F("hRawSumQZNC", "Raw sumQ ZNC;Amplitude [ADC counts];Counts",100, 0., 1200.);
204   TH1F * hRawSumQZNA = new TH1F("hRawSumQZNA", "Raw sumQ ZNA;Amplitude [ADC counts];Counts",100, 0., 1200.);
205   TH1F * hRawSumQZPC = new TH1F("hRawSumQZPC", "Raw sumQ ZPC;Amplitude [ADC counts];Counts",100, 0., 1200.);
206   TH1F * hRawSumQZPA = new TH1F("hRawSumQZPA", "Raw sumQ ZPA;Amplitude [ADC counts];Counts",100, 0., 1200.);
207   Add2RawsList(hRawSumQZNC, 14, expert, !image);
208   Add2RawsList(hRawSumQZNA, 15, expert, !image);
209   Add2RawsList(hRawSumQZPC, 16, expert, !image);
210   Add2RawsList(hRawSumQZPA, 17, expert, !image);
211   
212   TH1F * hRawTDCZEM1 = new TH1F("hRawTDCZEM1", "Raw TDC ZEM1;TDC [ns]",160, -350., -310.);
213   Add2RawsList(hRawTDCZEM1, 18, expert, !image);
214   TH1F * hRawTDCZPC = new TH1F("hRawTDCZPC", "Raw TDC ZPC;TDC [ns]",160, -350., -310.);
215   Add2RawsList(hRawTDCZPC, 19, expert, !image);
216   
217   TProfile * hRawADCProfs = new TProfile("hRawADCProfs", "ADC profiles;ADC id;Mean ADC values",22,-0.5,21.5,10.,1210.,"");
218   Add2RawsList(hRawADCProfs, 20, expert, !image);
219   TProfile * hRawTDCProfs = new TProfile("hRawTDCProfs", "TDC profiles;TDC id;Mean TDC values",6,0.5,6.5,-350.,-310.,"S");
220   Add2RawsList(hRawTDCProfs, 21, expert, !image);
221   
222   TH1F * hRawADCs = new TH1F("hRawADCs", "ADCs;ADC id;Mean ADC values",22,-0.5,21.5);
223   Add2RawsList(hRawADCs, 22, !expert, image);
224  
225   TH1F * hRawTDCs = new TH1F("hRawTDCs", "TDCs;TDC id;Mean TDC values",6,0.5,6.5);
226   hRawTDCs->SetMaximum(-310);
227   Add2RawsList(hRawTDCs, 23, !expert, image);
228   
229   TH2F *hZNCrawCentr  = new TH2F("hZNCrawCentr", "Centroid in ZNC;X (cm);Y(cm)", 100, -5.,5.,100,-5.,5.);
230   Add2RawsList(hZNCrawCentr, 24, expert, !image);
231   TH2F *hZNArawCentr  = new TH2F("hZNArawCentr", "Centroid in ZNA;X (cm);Y(cm)", 100, -5.,5.,100,-5.,5.);
232   Add2RawsList(hZNArawCentr, 25, expert, !image);
233   
234   TH2F *hTimeZDC = new TH2F("fHistTimeZDC", "ZDC;TDC timing A+C vs C-A; events", 120,-30,30,120,-600,-540);
235   Add2RawsList(hTimeZDC, !expert, image);
236 }
237
238 //____________________________________________________________________________
239 void AliZDCQADataMakerRec::InitRecPoints()
240 {
241   // create Digits histograms in Digits subdir
242   const Bool_t expert = kTRUE ; 
243   const Bool_t image  = kTRUE ; 
244   //
245   // ------------------- HIGH GAIN CHAIN ---------------------------
246   TH1F * hRecZNCTot = new TH1F("hRecZNCTot", "Rec signal in ZNC;Amplitude [ADC counts];Counts", 100, 0., 2000.);
247   TH1F * hRecZNATot = new TH1F("hRecZNATot", "Rec signal in ZNA;Amplitude [ADC counts];Counts", 100, 0., 2000.);
248   TH1F * hRecZPCTot = new TH1F("hRecZPCTot", "Rec signal in ZPC;Amplitude [ADC counts];Counts", 100, 0., 8000.);
249   TH1F * hRecZPATot = new TH1F("hRecZPATot", "Rec signal in ZPA;Amplitude [ADC counts];Counts", 100, 0., 8000.);
250   Add2RecPointsList(hRecZNCTot, 0, !expert, image);
251   Add2RecPointsList(hRecZNATot, 1, !expert, image);
252   Add2RecPointsList(hRecZPCTot, 2, !expert, image);
253   Add2RecPointsList(hRecZPATot, 3, !expert, image);
254   //
255   TH1F * hRecSumQZNC = new TH1F("hRecSumQZNC", "Rec summed 4 ZNC quadrants;Amplitude [ADC counts];Counts",100, 0., 2000.);
256   TH1F * hRecSumQZNA = new TH1F("hRecSumQZNA", "Rec summed 4 ZNA quadrants;Amplitude [ADC counts];Counts",100, 0., 2000.);
257   TH1F * hRecSumQZPC = new TH1F("hRecSumQZPC", "Rec summed 4 ZPC quadrants;Amplitude [ADC counts];Counts",100, 0., 2000.);
258   TH1F * hRecSumQZPA = new TH1F("hRecSumQZPA", "Rec summed 4 ZPA quadrants;Amplitude [ADC counts];Counts",100, 0., 2000.);
259   Add2RecPointsList(hRecSumQZNC, 4, expert, !image);
260   Add2RecPointsList(hRecSumQZNA, 5, expert, !image);
261   Add2RecPointsList(hRecSumQZPC, 6, expert, !image);
262   Add2RecPointsList(hRecSumQZPA, 7, expert, !image);
263   //
264   TH1F * hRecPMCZNC = new TH1F("hRecPMCZNC", "Rec common ZNC PMT;Amplitude [ADC counts];Counts",100, 0., 2000.);
265   TH1F * hRecPMCZNA = new TH1F("hRecPMCZNA", "Rec common ZNA PMT;Amplitude [ADC counts];Counts",100, 0., 2000.);
266   TH1F * hRecPMCZPC = new TH1F("hRecPMCZPC", "Rec common ZPC PMT;Amplitude [ADC counts];Counts",100, 0., 2000.);
267   TH1F * hRecPMCZPA = new TH1F("hRecPMCZPA", "Rec common ZPA PMT;Amplitude [ADC counts];Counts",100, 0., 2000.);
268   Add2RecPointsList(hRecPMCZNC, 8 , expert, !image);
269   Add2RecPointsList(hRecPMCZNA, 9 , expert, !image);
270   Add2RecPointsList(hRecPMCZPC, 10, expert, !image);
271   Add2RecPointsList(hRecPMCZPA, 11, expert, !image); 
272 }
273
274 //____________________________________________________________________________
275 void AliZDCQADataMakerRec::InitESDs()
276 {
277   //Booking ESDs histograms
278   //
279   const Bool_t expert = kTRUE ; 
280   const Bool_t image  = kTRUE ; 
281   
282   TH2F * hZNCcentr  = new TH2F("hZNCcentr", "Centroid in ZNC;X (cm);Y(cm)", 100, -5.,5.,100,-5.,5.);
283   TH2F * hZNAcentr  = new TH2F("hZNAcentr", "Centroid in ZNA;X (cm);Y(cm)", 100, -5.,5.,100,-5.,5.);
284   Add2ESDsList(hZNCcentr, 0, !expert, image);
285   Add2ESDsList(hZNAcentr, 1, !expert, image);
286   //
287   // ------------------- HIGH GAIN CHAIN ---------------------------
288   TH1F * hESDZNCTot = new TH1F("hESDZNCTot", "Energy in ZNC", 100, 0., 4000.);
289   TH1F * hESDZNATot = new TH1F("hESDZNATot", "Energy in ZNA", 100, 0., 4000.);
290   TH1F * hESDZPCTot = new TH1F("hESDZPCTot", "Energy in ZPC", 100, 0., 4000.);
291   TH1F * hESDZPATot = new TH1F("hESDZPATot", "Energy in ZPA", 100, 0., 4000.);
292   Add2ESDsList(hESDZNCTot, 2, !expert, image);
293   Add2ESDsList(hESDZNATot, 3, !expert, image);
294   Add2ESDsList(hESDZPCTot, 4, !expert, image);
295   Add2ESDsList(hESDZPATot, 5, !expert, image);
296   //
297   TH1F * hESDZEM1 = new TH1F("hESDZEM1", "Energy in ZEM1", 100, 0., 2000.);
298   TH1F * hESDZEM2 = new TH1F("hESDZEM2", "Energy in ZEM2", 100, 0., 2000.);
299   Add2ESDsList(hESDZEM1,6, !expert, image);
300   Add2ESDsList(hESDZEM2,7, !expert, image);
301   //
302   TH1F * hESDSumQZNC = new TH1F("hESDSumQZNC", "Sum of 4 ZNC energy",100, 0., 2000.);
303   TH1F * hESDSumQZNA = new TH1F("hESDSumQZNA", "Sum of 4 ZNA energy",100, 0., 2000.);
304   TH1F * hESDSumQZPC = new TH1F("hESDSumQZPC", "Sum of 4 ZPC energy",100, 0., 2000.);
305   TH1F * hESDSumQZPA = new TH1F("hESDSumQZPA", "Sum of 4 ZPA energy",100, 0., 2000.);
306   Add2ESDsList(hESDSumQZNC, 8, expert, !image);
307   Add2ESDsList(hESDSumQZNA, 9, expert, !image);
308   Add2ESDsList(hESDSumQZPC, 10, expert, !image);
309   Add2ESDsList(hESDSumQZPA, 11, expert, !image);
310   //
311   TH1F * hESDPMCZNC = new TH1F("hESDPMCZNC", "Energy in ZNC PMC",100, 0., 2000.);
312   TH1F * hESDPMCZNA = new TH1F("hESDPMCZNA", "Energy in ZNA PMC",100, 0., 2000.);
313   TH1F * hESDPMCZPC = new TH1F("hESDPMCZPC", "Energy in ZPC PMC",100, 0., 2000.);
314   TH1F * hESDPMCZPA = new TH1F("hESDPMCZPA", "Energy in ZPA PMC",100, 0., 2000.);
315   Add2ESDsList(hESDPMCZNC, 12, expert, !image);
316   Add2ESDsList(hESDPMCZNA, 13, expert, !image);
317   Add2ESDsList(hESDPMCZPC, 14, expert, !image);
318   Add2ESDsList(hESDPMCZPA, 15, expert, !image);
319   // 
320   // ------------------- LOW GAIN CHAIN ---------------------------
321   TH1F * hESDSumQZNClg = new TH1F("hESDSumQZNClg", "Sum of 4 lg ZNC sectors",100, 0., 4000.);
322   TH1F * hESDSumQZNAlg = new TH1F("hESDSumQZNAlg", "Sum of 4 lg ZNA sectors",100, 0., 4000.);
323   TH1F * hESDSumQZPClg = new TH1F("hESDSumQZPClg", "Sum of 4 lg ZPC sectors",100, 0., 4000.);
324   TH1F * hESDSumQZPAlg = new TH1F("hESDSumQZPAlg", "Sum of 4 lg ZPA sectors",100, 0., 4000.);
325   Add2ESDsList(hESDSumQZNClg, 16, expert, !image);
326   Add2ESDsList(hESDSumQZNAlg, 17, expert, !image);
327   Add2ESDsList(hESDSumQZPClg, 18, expert, !image);
328   Add2ESDsList(hESDSumQZPAlg, 19, expert, !image);
329   //
330   TH1F * hESDPMCZNClg = new TH1F("hESDPMCZNClg", "Signal in common ZNC lg PMT",100, 0., 4000.);
331   TH1F * hESDPMCZNAlg = new TH1F("hESDPMCZNAlg", "Signal in common ZNA lg PMT",100, 0., 4000.);
332   TH1F * hESDPMCZPClg = new TH1F("hESDPMCZPClg", "Signal in common ZPC lg PMT",100, 0., 4000.);
333   TH1F * hESDPMCZPAlg = new TH1F("hESDPMCZPAlg", "Signal in common ZPA lg PMT",100, 0., 4000.);
334   Add2ESDsList(hESDPMCZNClg, 20, expert, !image);
335   Add2ESDsList(hESDPMCZNAlg, 21, expert, !image);
336   Add2ESDsList(hESDPMCZPClg, 22, expert, !image);
337   Add2ESDsList(hESDPMCZPAlg, 23, expert, !image);
338 }
339
340 //___________________________________________________________________________
341 void AliZDCQADataMakerRec::MakeDigits(TTree *digitTree)
342 {
343   // makes data from Digit Tree
344   if(!GetDigitsData(0)) InitDigits();
345
346   if(!digitTree){
347     AliError("Can't get ZDC digit tree!!");
348     return; 
349   }     
350    
351   TBranch * branch = digitTree->GetBranch("ZDC");
352   if(!branch){
353     AliError("ZDC branch in digit tree not found"); 
354     return;
355   } 
356     
357   AliZDCDigit *digit = 0x0;
358   branch->SetAddress(&digit);
359      
360   Float_t adcSum_ZNC=0., adcSum_ZNA=0., adcSum_ZPC=0., adcSum_ZPA=0.;
361   Float_t adcSumQ_ZNC=0., adcSumQ_ZNA=0., adcSumQ_ZPC=0., adcSumQ_ZPA=0.;
362   Float_t adcSum_ZNC_lg=0., adcSum_ZNA_lg=0., adcSum_ZPC_lg=0., adcSum_ZPA_lg=0.;
363   Float_t adcSumQ_ZNC_lg=0., adcSumQ_ZNA_lg=0., adcSumQ_ZPC_lg=0., adcSumQ_ZPA_lg=0.;
364   
365   Int_t ndig = digitTree->GetEntries();
366   for(Int_t i=0; i<ndig; i++){
367       branch->GetEntry(i);
368       
369       if(digit->GetSector(0)==1 && digit->GetSector(1)!=5){
370           adcSum_ZNC += digit->GetADCValue(0);
371           adcSum_ZNC_lg += digit->GetADCValue(1);
372           //
373           if(digit->GetSector(1)!=0){
374               adcSumQ_ZNC += digit->GetADCValue(0);
375               adcSumQ_ZNC_lg+= digit->GetADCValue(1);
376           }
377           else{
378               GetDigitsData(8)->Fill(digit->GetADCValue(0));
379               GetDigitsData(20)->Fill(digit->GetADCValue(1));
380           }
381       }
382       else if(digit->GetSector(0)==2){
383           adcSum_ZPC += digit->GetADCValue(0);
384           adcSum_ZPC_lg += digit->GetADCValue(1);
385           //
386           if(digit->GetSector(1)!=0){
387               adcSumQ_ZPC += digit->GetADCValue(0);
388               adcSumQ_ZPC_lg+= digit->GetADCValue(1);
389           }
390           else{
391               GetDigitsData(10)->Fill(digit->GetADCValue(0));
392               GetDigitsData(22)->Fill(digit->GetADCValue(1));
393           }
394       }
395       else if(digit->GetSector(0)==4 && digit->GetSector(1)!=5){
396           adcSum_ZNA += digit->GetADCValue(0);
397           adcSum_ZNA_lg += digit->GetADCValue(1);
398           //
399           if(digit->GetSector(1)!=0){
400               adcSumQ_ZNA += digit->GetADCValue(0);
401               adcSumQ_ZNA_lg+= digit->GetADCValue(1);
402           }
403           else{
404               GetDigitsData(9)->Fill(digit->GetADCValue(0));
405               GetDigitsData(21)->Fill(digit->GetADCValue(1));
406           }
407       }
408       else if(digit->GetSector(0)==5){
409           adcSum_ZPA += digit->GetADCValue(0);
410           adcSum_ZPA_lg += digit->GetADCValue(1);
411           //
412           if(digit->GetSector(1)!=0){
413               adcSumQ_ZPA += digit->GetADCValue(0);
414               adcSumQ_ZPA_lg+= digit->GetADCValue(1);
415           }
416           else{
417               GetDigitsData(11)->Fill(digit->GetADCValue(0));
418               GetDigitsData(23)->Fill(digit->GetADCValue(1));
419           }
420       }
421   }
422   //
423   GetDigitsData(0)->Fill(adcSum_ZNC);
424   GetDigitsData(1)->Fill(adcSum_ZNA);
425   GetDigitsData(2)->Fill(adcSum_ZPC);
426   GetDigitsData(3)->Fill(adcSum_ZPA);
427   //
428   GetDigitsData(4)->Fill(adcSumQ_ZNC);
429   GetDigitsData(5)->Fill(adcSumQ_ZNA);
430   GetDigitsData(6)->Fill(adcSumQ_ZPC);
431   GetDigitsData(7)->Fill(adcSumQ_ZPA);
432   //
433   GetDigitsData(12)->Fill(adcSum_ZNC_lg);
434   GetDigitsData(13)->Fill(adcSum_ZNA_lg);
435   GetDigitsData(14)->Fill(adcSum_ZPC_lg);
436   GetDigitsData(15)->Fill(adcSum_ZPA_lg);
437   //
438   GetDigitsData(16)->Fill(adcSumQ_ZNC_lg);
439   GetDigitsData(17)->Fill(adcSumQ_ZNA_lg);
440   GetDigitsData(18)->Fill(adcSumQ_ZPC_lg);
441   GetDigitsData(19)->Fill(adcSumQ_ZPA_lg);
442   
443   delete digit;
444   digit=0;
445
446 }
447
448
449 //____________________________________________________________________________
450 void AliZDCQADataMakerRec::MakeRaws(AliRawReader *rawReader)
451 {
452   // Filling Raws QA histos
453   //
454   // Checking the event type 
455 //  if (rawReader->GetType()!=7){
456   
457     // Check if histograms already created for this Event Specie
458     if(!GetRawsData(0)) InitRaws();
459   
460     // Parameters for mean value pedestal subtraction
461     int const kNch = 24;
462     Float_t meanPed[2*kNch];    
463     for(Int_t jj=0; jj<2*kNch; jj++) meanPed[jj] = fPedCalibData->GetMeanPed(jj);
464     
465     Float_t zncSignal=0., znaSignal=0., zpcSignal=0., zpaSignal=0.;
466     Float_t zncSumQ=0., znaSumQ=0., zpcSumQ=0., zpaSumQ=0.;
467     Float_t zncpmC=0., znapmC=0., zpcpmC=0., zpapmC=0.;
468     Bool_t isZNCFired=kFALSE, isZPCFired=kFALSE, isZNAFired=kFALSE, isZPAFired=kFALSE;
469     Int_t  indZNC=0, indZNA=0, indZPC=0, indZPA=0;
470     Float_t zncTDC[10], zpcTDC[10], zem1TDC[10], zem2TDC[10], znaTDC[10], zpaTDC[10];
471     Float_t zncSumTDC[10], znaSumTDC[10];
472     for(Int_t i=0; i<10; i++){
473        zncTDC[i]=zpcTDC[i]=zem1TDC[i]=zem2TDC[i]=znaTDC[i]=zpaTDC[i]=zncSumTDC[i]=znaSumTDC[i]=-999.;
474     }
475     Float_t tdcGate=-999.;
476     Int_t iMultZNCTDC=0, iMultZPCTDC=0, iMultZEM1TDC=0, iMultZEM2TDC=0, iMultZNATDC=0, iMultZPATDC=0;
477     Int_t iMultTDCC=0, iMultTDCA=0;
478     //
479     const Float_t refSum = -568.5;
480     const Float_t refDelta = -2.1;
481     const Float_t sigmaSum = 3.25;
482     const Float_t sigmaDelta = 2.25;
483     
484     const Float_t x[4] = {-1.75, 1.75, -1.75, 1.75};
485     const Float_t y[4] = {-1.75, -1.75, 1.75, 1.75};
486     const Float_t alpha=0.5;
487     Float_t numXZNC=0., numYZNC=0., denZNC=0., wZNC=0.; 
488     Float_t numXZNA=0., numYZNA=0., denZNA=0., wZNA=0.; 
489     
490     rawReader->Reset();
491     AliZDCRawStream stream(rawReader);
492     while(stream.Next()){
493
494       if(stream.IsADCDataWord() && 
495          (stream.GetADCModule()==0 || stream.GetADCModule()==1)){
496        
497          Int_t det = stream.GetSector(0);
498          Int_t quad = stream.GetSector(1);
499          Int_t gain = stream.GetADCGain();
500          Int_t pedindex=0;
501          
502          // Stuff for pedestal subtraction
503          if(quad != 5){ // ZDCs (not reference PTMs)
504           Float_t pedSubVal=-99.;
505           if(det == 1){    
506             pedindex = quad;
507             if(gain == 0){
508               pedSubVal = (Float_t) (stream.GetADCValue()-meanPed[pedindex]); 
509               zncSignal  += pedSubVal; 
510               isZNCFired = kTRUE;
511               if(quad!=0){
512                 zncSumQ += pedSubVal;
513                 if(pedSubVal>0.&& zncpmC>50.){
514                   wZNC = TMath::Power(pedSubVal, alpha);
515                   numXZNC += x[quad-1]*wZNC;
516                   numYZNC += y[quad-1]*wZNC;
517                   denZNC += wZNC;
518                 }
519               }
520               else{
521                 zncpmC = pedSubVal;
522                 GetRawsData(10)->Fill(zncpmC);
523               }
524               indZNC++;
525               
526               GetRawsData(20)->Fill(pedindex, pedSubVal);
527             }
528           }
529           else if(det == 2){ 
530             pedindex = quad+5;
531             if(gain == 0){
532               pedSubVal = (Float_t) (stream.GetADCValue()-meanPed[pedindex]); 
533               zpcSignal += pedSubVal; 
534               isZPCFired = kTRUE;
535               if(quad!=0) zpcSumQ += pedSubVal;
536               else{
537                 zpcpmC = pedSubVal;
538                 GetRawsData(12)->Fill(zpcpmC);
539               }
540               indZPC++;
541               
542               GetRawsData(20)->Fill(pedindex, pedSubVal);
543             }
544           }
545           else if(det == 3){ 
546             pedindex = quad+9;
547             if(quad==1){     
548               if(gain == 0){
549                 pedSubVal = (Float_t) (stream.GetADCValue()-meanPed[pedindex]); 
550                 GetRawsData(4)->Fill(pedSubVal);
551               
552                 GetRawsData(20)->Fill(pedindex, pedSubVal);
553               }
554             }
555             else if(quad==2){ 
556               if(gain == 0){
557                 pedSubVal = (Float_t) (stream.GetADCValue()-meanPed[pedindex]); 
558                 GetRawsData(5)->Fill(pedSubVal); 
559               
560                 GetRawsData(20)->Fill(pedindex, pedSubVal);
561               }
562             }
563           }
564           else if(det == 4){       
565             pedindex = quad+12;
566             if(gain == 0){
567               pedSubVal = (Float_t) (stream.GetADCValue()-meanPed[pedindex]); 
568               znaSignal  += pedSubVal; 
569               isZNAFired = kTRUE;
570               if(quad!=0){
571                 znaSumQ += pedSubVal;
572                 if(pedSubVal>0.&& znapmC>50.) {
573                   wZNA = TMath::Power(pedSubVal, alpha);
574                   numXZNA += x[quad-1]*wZNA;
575                   numYZNA += y[quad-1]*wZNA;
576                   denZNA += wZNA;
577                 }
578               }
579               else{
580                 znapmC = pedSubVal;
581                 GetRawsData(11)->Fill(znapmC);
582               }
583               indZNA++;
584               
585               GetRawsData(20)->Fill(pedindex, pedSubVal);
586               //GetRawsData(22)->SetBinContent(pedindex+1, GetRawsData(20)->GetBinContent(pedindex+1));
587               //GetRawsData(22)->SetBinError(pedindex+1, GetRawsData(20)->GetBinError(pedindex+1));
588             }
589           }
590           else if(det == 5){
591             pedindex = quad+17;
592             if(gain == 0){
593               pedSubVal = (Float_t) (stream.GetADCValue()-meanPed[pedindex]); 
594               zpaSignal  += pedSubVal; 
595               isZPAFired = kTRUE;
596               if(quad!=0) zpaSumQ += pedSubVal;
597               else{
598                 zpapmC = pedSubVal;
599                 GetRawsData(13)->Fill(zpapmC);
600               }
601               indZPA++;
602               
603               GetRawsData(20)->Fill(pedindex, pedSubVal);
604               //GetRawsData(22)->SetBinContent(pedindex+1, GetRawsData(20)->GetBinContent(pedindex+1));
605               //GetRawsData(22)->SetBinError(pedindex+1, GetRawsData(20)->GetBinError(pedindex+1));
606             }
607           }
608                          
609          }
610
611          if(isZNCFired && indZNC==5){
612            GetRawsData(0)->Fill(zncSignal);
613            GetRawsData(6)->Fill(zncpmC, zncSumQ);
614            GetRawsData(14)->Fill(zncSumQ); 
615            //
616            Float_t xZNC, yZNC;           
617            if(denZNC!=0){
618              xZNC = numXZNC/denZNC;
619              yZNC = numYZNC/denZNC;
620            } 
621            else xZNC = yZNC = 999.;
622            GetRawsData(24)->Fill(xZNC, yZNC);
623          }
624          if(isZPCFired && indZPC==5){
625            GetRawsData(2)->Fill(zpcSignal);
626            GetRawsData(8)->Fill(zpcpmC, zpcSumQ);
627            GetRawsData(16)->Fill(zpcSumQ); 
628          }
629          if(isZNAFired && indZNA==5){ 
630            GetRawsData(1)->Fill(znaSignal);
631            GetRawsData(7)->Fill(znapmC, znaSumQ);
632            GetRawsData(15)->Fill(znaSumQ); 
633            //
634            Float_t xZNA, yZNA;
635            if(denZNA!=0){
636              xZNA = numXZNA/denZNA;
637              yZNA = numYZNA/denZNA;
638            } 
639            else xZNA = yZNA = 999.;
640            GetRawsData(25)->Fill(xZNA, yZNA);
641          }
642          if(isZPAFired && indZPA==5){ 
643            GetRawsData(3)->Fill(zpaSignal);
644            GetRawsData(9)->Fill(zpapmC, zpaSumQ);
645            GetRawsData(17)->Fill(zpaSumQ); 
646          }
647          
648          if(indZNC==5){
649            zncSignal = zncSumQ = zncpmC = 0;
650            isZNCFired=kFALSE; indZNC=0;
651          }
652          if(indZPC==5){
653            zpcSignal = zpcSumQ = zpcpmC = 0;
654            isZPCFired=kFALSE; indZPC=0;
655          }
656          if(indZNA==5){
657            znaSignal = znaSumQ = znapmC = 0;
658            isZNAFired=kFALSE; indZNA=0;
659          }
660          if(indZPA==5){
661            zpaSignal = zpaSumQ = zpapmC = 0;
662            isZPAFired=kFALSE; indZPA=0;
663          }
664          
665       } //IsADCDataWord && signal ADCs
666       else if(stream.IsZDCTDCDatum()){
667          if(stream.GetChannel()==1){
668             zncTDC[iMultZNCTDC] = (stream.GetZDCTDCDatum());
669             iMultZNCTDC++;
670          }
671          else if(stream.GetChannel()==3){
672             zpcTDC[iMultZPCTDC] = (stream.GetZDCTDCDatum());
673             iMultZPCTDC++;
674          }
675          else if(stream.GetChannel()==5){
676             znaTDC[iMultZNATDC] = (stream.GetZDCTDCDatum());
677             iMultZNATDC++;
678          }
679          else if(stream.GetChannel()==7){
680             zpaTDC[iMultZPATDC] = (stream.GetZDCTDCDatum());
681             iMultZPATDC++;
682          }
683          else if(stream.GetChannel()==8){
684             zem1TDC[iMultZEM1TDC] = (stream.GetZDCTDCDatum());
685             iMultZEM1TDC++;
686          }
687          else if(stream.GetChannel()==9){
688             zem2TDC[iMultZEM2TDC] = (stream.GetZDCTDCDatum());
689             iMultZEM2TDC++;
690          }
691          else if(stream.GetChannel()==10){
692             zncSumTDC[iMultZEM2TDC] = (stream.GetZDCTDCDatum());
693             iMultTDCC++;
694          }
695          else if(stream.GetChannel()==12){
696             znaSumTDC[iMultZEM2TDC] = (stream.GetZDCTDCDatum());
697             iMultTDCA++;
698          }
699          else if(stream.GetChannel()==14) tdcGate = (stream.GetZDCTDCDatum());
700          
701          if(stream.GetChannel()==16 && tdcGate!=-999.){
702            for(Int_t iHit=0; iHit<10; iHit++){
703               if(zncTDC[iHit]!=-999.){
704                 GetRawsData(21)->Fill(1, zncTDC[iHit]-tdcGate);
705                 //GetRawsData(23)->SetBinContent(1, GetRawsData(21)->GetBinContent(1));
706               }
707               if(zpcTDC[iHit]!=-999.){
708                 Float_t diffZPC = zpcTDC[iHit]-tdcGate;
709                 GetRawsData(19)->Fill(diffZPC);
710                 GetRawsData(21)->Fill(2, diffZPC);
711                 //GetRawsData(23)->SetBinContent(2,  GetRawsData(21)->GetBinContent(3));
712               }
713               if(znaTDC[iHit]!=-999.){
714                 GetRawsData(21)->Fill(3, znaTDC[iHit]-tdcGate);
715                 //GetRawsData(23)->SetBinContent(3,  GetRawsData(21)->GetBinContent(5));
716               }
717               if(zpaTDC[iHit]!=-999.){
718                 GetRawsData(21)->Fill(4, zpaTDC[iHit]-tdcGate);
719                 //GetRawsData(23)->SetBinContent(4,  GetRawsData(21)->GetBinContent(7));
720               }
721               if(zem1TDC[iHit]!=-999.){
722                 Float_t diffZEM1 = zem1TDC[iHit]-tdcGate;
723                 GetRawsData(18)->Fill(diffZEM1);
724                 GetRawsData(21)->Fill(5, diffZEM1);
725                 //GetRawsData(23)->SetBinContent(5,  GetRawsData(21)->GetBinContent(8));
726               }
727               if(zem2TDC[iHit]!=-999.){
728                 GetRawsData(21)->Fill(6, zem2TDC[iHit]-tdcGate);
729                 //GetRawsData(23)->SetBinContent(6,  GetRawsData(21)->GetBinContent(9));
730               }
731               if(zncSumTDC[iHit]!=-999.){
732                  Float_t tdcC = zncSumTDC[iHit]-tdcGate;
733                  if(znaSumTDC[iHit]!=-999.){
734                     Float_t tdcA = znaSumTDC[iHit]-tdcGate;
735                     if (((tdcC-tdcA-refDelta)*(tdcC-tdcA-refDelta)/(sigmaDelta*sigmaDelta) +
736                         (tdcC+tdcA-refSum)*(tdcC+tdcA-refSum)/(sigmaSum*sigmaSum))< 1.0)
737                         GetRawsData(26)->Fill(tdcC-tdcA,tdcC+tdcA);
738                     
739                  }
740               }
741            }
742            //
743            tdcGate = -999.;
744            for(Int_t i=0; i<10; i++){
745               zpcTDC[i] = zem1TDC[i] = zncSumTDC[i] = znaSumTDC[i] = -999.;
746            } 
747          }
748       }
749     
750     } //stream.Next()
751 //  } // check on event type
752 //  else{
753 //    AliDebug(1,Form("Skipping non-physics event for QA -> event type %d \n", rawReader->GetType())); 
754 //  }
755 }
756
757 //____________________________________________________________________________
758 void AliZDCQADataMakerRec::MakeRecPoints(TTree * clustersTree)
759 {
760   // Filling QA histos from RecPoints
761
762   TBranch *branch = clustersTree->GetBranch("ZDC");
763   if(!branch){ 
764     AliError("Can't get the ZDC branch for rec points!");
765     return;
766   }
767   
768   if(!GetRecPointsData(0)) InitRecPoints() ;
769
770   Float_t sum_ZNC=0., sum_ZNA=0., sum_ZPC=0., sum_ZPA=0.;
771   Float_t sumQ_ZNC=0., sumQ_ZNA=0., sumQ_ZPC=0., sumQ_ZPA=0.;
772
773   AliZDCReco reco;
774   AliZDCReco* preco = &reco;
775   clustersTree->SetBranchAddress("ZDC", &preco);
776
777   clustersTree->GetEntry(0);
778   for(Int_t i=0; i<5; i++){
779     sum_ZNC += reco.GetZN1HREnTow(i);
780     sum_ZPC += reco.GetZN2HREnTow(i);
781     sum_ZNA += reco.GetZP1HREnTow(i);
782     sum_ZPA += reco.GetZP2HREnTow(i);
783     if(i==0){
784       GetRecPointsData(8)->Fill(reco.GetZN1HREnTow(i));
785       GetRecPointsData(9)->Fill(reco.GetZN2HREnTow(i));
786       GetRecPointsData(10)->Fill(reco.GetZP1HREnTow(i));
787       GetRecPointsData(11)->Fill(reco.GetZP2HREnTow(i));
788     }
789     else{
790       sumQ_ZNC += reco.GetZN1HREnTow(i);
791       sumQ_ZPC += reco.GetZN2HREnTow(i);
792       sumQ_ZNA += reco.GetZP1HREnTow(i);
793       sumQ_ZPA += reco.GetZP2HREnTow(i);
794     }
795   }
796   
797   GetRecPointsData(0)->Fill(sum_ZNC);
798   GetRecPointsData(1)->Fill(sum_ZNA);
799   GetRecPointsData(2)->Fill(sum_ZPC);
800   GetRecPointsData(3)->Fill(sum_ZPA);
801   //
802   GetRecPointsData(4)->Fill(sumQ_ZNC);
803   GetRecPointsData(5)->Fill(sumQ_ZNA);
804   GetRecPointsData(6)->Fill(sumQ_ZPC);
805   GetRecPointsData(7)->Fill(sumQ_ZPA);
806   
807 }  
808
809 //____________________________________________________________________________
810 void AliZDCQADataMakerRec::MakeESDs(AliESDEvent * esd)
811 {
812   // make QA data from ESDs
813   //
814   
815   // Check id histograms already created for this Event Specie
816   if(!GetESDsData(0)) InitESDs() ;
817
818   AliESDZDC * zdcESD =  esd->GetESDZDC();
819   //
820   TString beamType = esd->GetBeamType();
821   Double_t centr_ZNC[2]={999.,999}, centr_ZNA[2]={999.,999};
822   if(((beamType.CompareTo("pp"))==0) || ((beamType.CompareTo("p-p"))==0)
823      ||((beamType.CompareTo("PP"))==0) || ((beamType.CompareTo("P-P"))==0)){
824     zdcESD->GetZNCentroidInpp(centr_ZNC, centr_ZNA);
825   }
826   else if((beamType.CompareTo("A-A")) == 0 || (beamType.CompareTo("Pb-Pb")) == 0){
827     Float_t beamEne = esd->GetBeamEnergy();
828     zdcESD->GetZNCentroidInPbPb(beamEne, centr_ZNC, centr_ZNA);
829   }
830   else printf("\n WARNING!!! AliZDCQADataMakerRec::MakeESDs: can't calculate centroids for beam type: %s\n\n",beamType.Data());
831   GetESDsData(0)->Fill(centr_ZNC[0], centr_ZNC[1]);
832   GetESDsData(1)->Fill(centr_ZNA[0], centr_ZNA[1]);
833
834   //
835   GetESDsData(2)->Fill(esd->GetZDCN1Energy());
836   GetESDsData(3)->Fill(esd->GetZDCN2Energy());
837   GetESDsData(4)->Fill(esd->GetZDCP1Energy());
838   GetESDsData(5)->Fill(esd->GetZDCP2Energy());
839   GetESDsData(6)->Fill(esd->GetZDCEMEnergy(0));
840   GetESDsData(7)->Fill(esd->GetZDCEMEnergy(1));
841   //
842   Double_t sumQZNC=0., sumQZPC=0., sumQZNA=0., sumQZPA=0.;
843   Double_t sumQZNC_lg=0., sumQZPC_lg=0., sumQZNA_lg=0., sumQZPA_lg=0.;
844   //
845   const Double_t *towZNC, *towZPC, *towZNA, *towZPA;
846   const Double_t *towZNC_lg, *towZPC_lg, *towZNA_lg, *towZPA_lg;
847   //
848   towZNC = zdcESD->GetZN1TowerEnergy();
849   towZPC = zdcESD->GetZP1TowerEnergy();
850   towZNA = zdcESD->GetZN2TowerEnergy();
851   towZPA = zdcESD->GetZP2TowerEnergy();
852   //
853   towZNC_lg = zdcESD->GetZN1TowerEnergyLR();
854   towZPC_lg = zdcESD->GetZP1TowerEnergyLR();
855   towZNA_lg = zdcESD->GetZN2TowerEnergyLR();
856   towZPA_lg = zdcESD->GetZP2TowerEnergyLR();
857   //
858   for(Int_t i=0; i<5; i++){
859      if(i==0){
860        GetESDsData(12)->Fill(towZNC[i]);
861        GetESDsData(13)->Fill(towZNA[i]);
862        GetESDsData(14)->Fill(towZPC[i]);
863        GetESDsData(15)->Fill(towZPA[i]);
864        //
865        GetESDsData(20)->Fill(towZNC_lg[i]);
866        GetESDsData(21)->Fill(towZNA_lg[i]);
867        GetESDsData(22)->Fill(towZPC_lg[i]);
868        GetESDsData(23)->Fill(towZPA_lg[i]);
869      }
870      else{
871        sumQZNC += towZNC[i];
872        sumQZPC += towZPC[i];
873        sumQZNA += towZNA[i];
874        sumQZPA += towZPA[i];
875        //
876        sumQZNC_lg += towZNC_lg[i];
877        sumQZPC_lg += towZPC_lg[i];
878        sumQZNA_lg += towZNA_lg[i];
879        sumQZPA_lg += towZPA_lg[i];
880      }
881   }
882   GetESDsData(8)->Fill(sumQZNC);
883   GetESDsData(9)->Fill(sumQZNA);
884   GetESDsData(10)->Fill(sumQZPC);
885   GetESDsData(11)->Fill(sumQZPA);
886   //
887   GetESDsData(16)->Fill(sumQZNC_lg);
888   GetESDsData(17)->Fill(sumQZNA_lg);
889   GetESDsData(18)->Fill(sumQZPC_lg);
890   GetESDsData(19)->Fill(sumQZPA_lg);
891 }
892
893 //____________________________________________________________________________
894 void AliZDCQADataMakerRec::StartOfDetectorCycle()
895 {
896   //Detector specific actions at start of cycle
897
898   fPedCalibData = GetPedCalibData();
899   
900 }
901
902 //____________________________________________________________________________ 
903 void AliZDCQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
904 {
905   //Detector specific actions at end of cycle
906   // do the QA checking
907   if( task == AliQAv1::kRAWS){
908      if (!GetRawsData(20) || !GetRawsData(21) || !GetRawsData(24)  || !GetRawsData(25)) {
909          printf("  WARNING!!! AliZDCQADataMaker Rec -> Not all histogram for DQM found!\n"); 
910      }
911      else{
912        for(Int_t ibin=1; ibin<=GetRawsData(20)->GetNbinsX(); ibin++){
913           GetRawsData(22)->SetBinContent(ibin, GetRawsData(20)->GetBinContent(ibin)); 
914           GetRawsData(22)->SetBinError(ibin, GetRawsData(20)->GetBinError(ibin));
915           GetRawsData(22)->SetLineColor(kBlue); GetRawsData(22)->SetLineWidth(2);
916        }
917        for(Int_t ibin=1; ibin<=GetRawsData(21)->GetNbinsX(); ibin++){
918           GetRawsData(23)->SetBinContent(ibin, GetRawsData(21)->GetBinContent(ibin)); 
919           GetRawsData(23)->SetBinError(ibin, GetRawsData(21)->GetBinError(ibin));
920           GetRawsData(23)->SetLineColor(kAzure-3); GetRawsData(23)->SetLineWidth(2);
921        }
922      }
923   }
924         
925   AliQAChecker::Instance()->Run(AliQAv1::kZDC, task, list) ;  
926 }