]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/AliZDCQADataMakerRec.cxx
Updated DAs
[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
235 //____________________________________________________________________________
236 void AliZDCQADataMakerRec::InitRecPoints()
237 {
238   // create Digits histograms in Digits subdir
239   const Bool_t expert = kTRUE ; 
240   const Bool_t image  = kTRUE ; 
241   //
242   // ------------------- HIGH GAIN CHAIN ---------------------------
243   TH1F * hRecZNCTot = new TH1F("hRecZNCTot", "Rec signal in ZNC;Amplitude [ADC counts];Counts", 100, 0., 2000.);
244   TH1F * hRecZNATot = new TH1F("hRecZNATot", "Rec signal in ZNA;Amplitude [ADC counts];Counts", 100, 0., 2000.);
245   TH1F * hRecZPCTot = new TH1F("hRecZPCTot", "Rec signal in ZPC;Amplitude [ADC counts];Counts", 100, 0., 8000.);
246   TH1F * hRecZPATot = new TH1F("hRecZPATot", "Rec signal in ZPA;Amplitude [ADC counts];Counts", 100, 0., 8000.);
247   Add2RecPointsList(hRecZNCTot, 0, !expert, image);
248   Add2RecPointsList(hRecZNATot, 1, !expert, image);
249   Add2RecPointsList(hRecZPCTot, 2, !expert, image);
250   Add2RecPointsList(hRecZPATot, 3, !expert, image);
251   //
252   TH1F * hRecSumQZNC = new TH1F("hRecSumQZNC", "Rec summed 4 ZNC quadrants;Amplitude [ADC counts];Counts",100, 0., 2000.);
253   TH1F * hRecSumQZNA = new TH1F("hRecSumQZNA", "Rec summed 4 ZNA quadrants;Amplitude [ADC counts];Counts",100, 0., 2000.);
254   TH1F * hRecSumQZPC = new TH1F("hRecSumQZPC", "Rec summed 4 ZPC quadrants;Amplitude [ADC counts];Counts",100, 0., 2000.);
255   TH1F * hRecSumQZPA = new TH1F("hRecSumQZPA", "Rec summed 4 ZPA quadrants;Amplitude [ADC counts];Counts",100, 0., 2000.);
256   Add2RecPointsList(hRecSumQZNC, 4, expert, !image);
257   Add2RecPointsList(hRecSumQZNA, 5, expert, !image);
258   Add2RecPointsList(hRecSumQZPC, 6, expert, !image);
259   Add2RecPointsList(hRecSumQZPA, 7, expert, !image);
260   //
261   TH1F * hRecPMCZNC = new TH1F("hRecPMCZNC", "Rec common ZNC PMT;Amplitude [ADC counts];Counts",100, 0., 2000.);
262   TH1F * hRecPMCZNA = new TH1F("hRecPMCZNA", "Rec common ZNA PMT;Amplitude [ADC counts];Counts",100, 0., 2000.);
263   TH1F * hRecPMCZPC = new TH1F("hRecPMCZPC", "Rec common ZPC PMT;Amplitude [ADC counts];Counts",100, 0., 2000.);
264   TH1F * hRecPMCZPA = new TH1F("hRecPMCZPA", "Rec common ZPA PMT;Amplitude [ADC counts];Counts",100, 0., 2000.);
265   Add2RecPointsList(hRecPMCZNC, 8 , expert, !image);
266   Add2RecPointsList(hRecPMCZNA, 9 , expert, !image);
267   Add2RecPointsList(hRecPMCZPC, 10, expert, !image);
268   Add2RecPointsList(hRecPMCZPA, 11, expert, !image); 
269 }
270
271 //____________________________________________________________________________
272 void AliZDCQADataMakerRec::InitESDs()
273 {
274   //Booking ESDs histograms
275   //
276   const Bool_t expert = kTRUE ; 
277   const Bool_t image  = kTRUE ; 
278   
279   TH2F * hZNCcentr  = new TH2F("hZNCcentr", "Centroid in ZNC;X (cm);Y(cm)", 100, -5.,5.,100,-5.,5.);
280   TH2F * hZNAcentr  = new TH2F("hZNAcentr", "Centroid in ZNA;X (cm);Y(cm)", 100, -5.,5.,100,-5.,5.);
281   Add2ESDsList(hZNCcentr, 0, !expert, image);
282   Add2ESDsList(hZNAcentr, 1, !expert, image);
283   //
284   // ------------------- HIGH GAIN CHAIN ---------------------------
285   TH1F * hESDZNCTot = new TH1F("hESDZNCTot", "Energy in ZNC", 100, 0., 4000.);
286   TH1F * hESDZNATot = new TH1F("hESDZNATot", "Energy in ZNA", 100, 0., 4000.);
287   TH1F * hESDZPCTot = new TH1F("hESDZPCTot", "Energy in ZPC", 100, 0., 4000.);
288   TH1F * hESDZPATot = new TH1F("hESDZPATot", "Energy in ZPA", 100, 0., 4000.);
289   Add2ESDsList(hESDZNCTot, 2, !expert, image);
290   Add2ESDsList(hESDZNATot, 3, !expert, image);
291   Add2ESDsList(hESDZPCTot, 4, !expert, image);
292   Add2ESDsList(hESDZPATot, 5, !expert, image);
293   //
294   TH1F * hESDZEM1 = new TH1F("hESDZEM1", "Energy in ZEM1", 100, 0., 2000.);
295   TH1F * hESDZEM2 = new TH1F("hESDZEM2", "Energy in ZEM2", 100, 0., 2000.);
296   Add2ESDsList(hESDZEM1,6, !expert, image);
297   Add2ESDsList(hESDZEM2,7, !expert, image);
298   //
299   TH1F * hESDSumQZNC = new TH1F("hESDSumQZNC", "Sum of 4 ZNC energy",100, 0., 2000.);
300   TH1F * hESDSumQZNA = new TH1F("hESDSumQZNA", "Sum of 4 ZNA energy",100, 0., 2000.);
301   TH1F * hESDSumQZPC = new TH1F("hESDSumQZPC", "Sum of 4 ZPC energy",100, 0., 2000.);
302   TH1F * hESDSumQZPA = new TH1F("hESDSumQZPA", "Sum of 4 ZPA energy",100, 0., 2000.);
303   Add2ESDsList(hESDSumQZNC, 8, expert, !image);
304   Add2ESDsList(hESDSumQZNA, 9, expert, !image);
305   Add2ESDsList(hESDSumQZPC, 10, expert, !image);
306   Add2ESDsList(hESDSumQZPA, 11, expert, !image);
307   //
308   TH1F * hESDPMCZNC = new TH1F("hESDPMCZNC", "Energy in ZNC PMC",100, 0., 2000.);
309   TH1F * hESDPMCZNA = new TH1F("hESDPMCZNA", "Energy in ZNA PMC",100, 0., 2000.);
310   TH1F * hESDPMCZPC = new TH1F("hESDPMCZPC", "Energy in ZPC PMC",100, 0., 2000.);
311   TH1F * hESDPMCZPA = new TH1F("hESDPMCZPA", "Energy in ZPA PMC",100, 0., 2000.);
312   Add2ESDsList(hESDPMCZNC, 12, expert, !image);
313   Add2ESDsList(hESDPMCZNA, 13, expert, !image);
314   Add2ESDsList(hESDPMCZPC, 14, expert, !image);
315   Add2ESDsList(hESDPMCZPA, 15, expert, !image);
316   // 
317   // ------------------- LOW GAIN CHAIN ---------------------------
318   TH1F * hESDSumQZNClg = new TH1F("hESDSumQZNClg", "Sum of 4 lg ZNC sectors",100, 0., 4000.);
319   TH1F * hESDSumQZNAlg = new TH1F("hESDSumQZNAlg", "Sum of 4 lg ZNA sectors",100, 0., 4000.);
320   TH1F * hESDSumQZPClg = new TH1F("hESDSumQZPClg", "Sum of 4 lg ZPC sectors",100, 0., 4000.);
321   TH1F * hESDSumQZPAlg = new TH1F("hESDSumQZPAlg", "Sum of 4 lg ZPA sectors",100, 0., 4000.);
322   Add2ESDsList(hESDSumQZNClg, 16, expert, !image);
323   Add2ESDsList(hESDSumQZNAlg, 17, expert, !image);
324   Add2ESDsList(hESDSumQZPClg, 18, expert, !image);
325   Add2ESDsList(hESDSumQZPAlg, 19, expert, !image);
326   //
327   TH1F * hESDPMCZNClg = new TH1F("hESDPMCZNClg", "Signal in common ZNC lg PMT",100, 0., 4000.);
328   TH1F * hESDPMCZNAlg = new TH1F("hESDPMCZNAlg", "Signal in common ZNA lg PMT",100, 0., 4000.);
329   TH1F * hESDPMCZPClg = new TH1F("hESDPMCZPClg", "Signal in common ZPC lg PMT",100, 0., 4000.);
330   TH1F * hESDPMCZPAlg = new TH1F("hESDPMCZPAlg", "Signal in common ZPA lg PMT",100, 0., 4000.);
331   Add2ESDsList(hESDPMCZNClg, 20, expert, !image);
332   Add2ESDsList(hESDPMCZNAlg, 21, expert, !image);
333   Add2ESDsList(hESDPMCZPClg, 22, expert, !image);
334   Add2ESDsList(hESDPMCZPAlg, 23, expert, !image);
335 }
336
337 //___________________________________________________________________________
338 void AliZDCQADataMakerRec::MakeDigits(TTree *digitTree)
339 {
340   // makes data from Digit Tree
341   if(!GetDigitsData(0)) InitDigits();
342
343   if(!digitTree){
344     AliError("Can't get ZDC digit tree!!");
345     return; 
346   }     
347    
348   TBranch * branch = digitTree->GetBranch("ZDC");
349   if(!branch){
350     AliError("ZDC branch in digit tree not found"); 
351     return;
352   } 
353     
354   AliZDCDigit *digit = 0x0;
355   branch->SetAddress(&digit);
356      
357   Float_t adcSum_ZNC=0., adcSum_ZNA=0., adcSum_ZPC=0., adcSum_ZPA=0.;
358   Float_t adcSumQ_ZNC=0., adcSumQ_ZNA=0., adcSumQ_ZPC=0., adcSumQ_ZPA=0.;
359   Float_t adcSum_ZNC_lg=0., adcSum_ZNA_lg=0., adcSum_ZPC_lg=0., adcSum_ZPA_lg=0.;
360   Float_t adcSumQ_ZNC_lg=0., adcSumQ_ZNA_lg=0., adcSumQ_ZPC_lg=0., adcSumQ_ZPA_lg=0.;
361   
362   Int_t ndig = digitTree->GetEntries();
363   for(Int_t i=0; i<ndig; i++){
364       branch->GetEntry(i);
365       
366       if(digit->GetSector(0)==1 && digit->GetSector(1)!=5){
367           adcSum_ZNC += digit->GetADCValue(0);
368           adcSum_ZNC_lg += digit->GetADCValue(1);
369           //
370           if(digit->GetSector(1)!=0){
371               adcSumQ_ZNC += digit->GetADCValue(0);
372               adcSumQ_ZNC_lg+= digit->GetADCValue(1);
373           }
374           else{
375               GetDigitsData(8)->Fill(digit->GetADCValue(0));
376               GetDigitsData(20)->Fill(digit->GetADCValue(1));
377           }
378       }
379       else if(digit->GetSector(0)==2){
380           adcSum_ZPC += digit->GetADCValue(0);
381           adcSum_ZPC_lg += digit->GetADCValue(1);
382           //
383           if(digit->GetSector(1)!=0){
384               adcSumQ_ZPC += digit->GetADCValue(0);
385               adcSumQ_ZPC_lg+= digit->GetADCValue(1);
386           }
387           else{
388               GetDigitsData(10)->Fill(digit->GetADCValue(0));
389               GetDigitsData(22)->Fill(digit->GetADCValue(1));
390           }
391       }
392       else if(digit->GetSector(0)==4 && digit->GetSector(1)!=5){
393           adcSum_ZNA += digit->GetADCValue(0);
394           adcSum_ZNA_lg += digit->GetADCValue(1);
395           //
396           if(digit->GetSector(1)!=0){
397               adcSumQ_ZNA += digit->GetADCValue(0);
398               adcSumQ_ZNA_lg+= digit->GetADCValue(1);
399           }
400           else{
401               GetDigitsData(9)->Fill(digit->GetADCValue(0));
402               GetDigitsData(21)->Fill(digit->GetADCValue(1));
403           }
404       }
405       else if(digit->GetSector(0)==5){
406           adcSum_ZPA += digit->GetADCValue(0);
407           adcSum_ZPA_lg += digit->GetADCValue(1);
408           //
409           if(digit->GetSector(1)!=0){
410               adcSumQ_ZPA += digit->GetADCValue(0);
411               adcSumQ_ZPA_lg+= digit->GetADCValue(1);
412           }
413           else{
414               GetDigitsData(11)->Fill(digit->GetADCValue(0));
415               GetDigitsData(23)->Fill(digit->GetADCValue(1));
416           }
417       }
418   }
419   //
420   GetDigitsData(0)->Fill(adcSum_ZNC);
421   GetDigitsData(1)->Fill(adcSum_ZNA);
422   GetDigitsData(2)->Fill(adcSum_ZPC);
423   GetDigitsData(3)->Fill(adcSum_ZPA);
424   //
425   GetDigitsData(4)->Fill(adcSumQ_ZNC);
426   GetDigitsData(5)->Fill(adcSumQ_ZNA);
427   GetDigitsData(6)->Fill(adcSumQ_ZPC);
428   GetDigitsData(7)->Fill(adcSumQ_ZPA);
429   //
430   GetDigitsData(12)->Fill(adcSum_ZNC_lg);
431   GetDigitsData(13)->Fill(adcSum_ZNA_lg);
432   GetDigitsData(14)->Fill(adcSum_ZPC_lg);
433   GetDigitsData(15)->Fill(adcSum_ZPA_lg);
434   //
435   GetDigitsData(16)->Fill(adcSumQ_ZNC_lg);
436   GetDigitsData(17)->Fill(adcSumQ_ZNA_lg);
437   GetDigitsData(18)->Fill(adcSumQ_ZPC_lg);
438   GetDigitsData(19)->Fill(adcSumQ_ZPA_lg);
439   
440   delete digit;
441   digit=0;
442
443 }
444
445
446 //____________________________________________________________________________
447 void AliZDCQADataMakerRec::MakeRaws(AliRawReader *rawReader)
448 {
449   // Filling Raws QA histos
450   //
451   // Checking the event type 
452 //  if (rawReader->GetType()!=7){
453   
454     // Check if histograms already created for this Event Specie
455     if(!GetRawsData(0)) InitRaws();
456   
457     // Parameters for mean value pedestal subtraction
458     int const kNch = 24;
459     Float_t meanPed[2*kNch];    
460     for(Int_t jj=0; jj<2*kNch; jj++) meanPed[jj] = fPedCalibData->GetMeanPed(jj);
461     
462     Float_t zncSignal=0., znaSignal=0., zpcSignal=0., zpaSignal=0.;
463     Float_t zncSumQ=0., znaSumQ=0., zpcSumQ=0., zpaSumQ=0.;
464     Float_t zncpmC=0., znapmC=0., zpcpmC=0., zpapmC=0.;
465     Bool_t isZNCFired=kFALSE, isZPCFired=kFALSE, isZNAFired=kFALSE, isZPAFired=kFALSE;
466     Int_t  indZNC=0, indZNA=0, indZPC=0, indZPA=0;
467     Float_t zncTDC[10], zpcTDC[10], zem1TDC[10], zem2TDC[10], znaTDC[10], zpaTDC[10];
468     for(Int_t i=0; i<10; i++){
469        zncTDC[i]=zpcTDC[i]=zem1TDC[i]=zem2TDC[i]=znaTDC[i]=zpaTDC[i]=-999.;
470     }
471     Float_t tdcGate=-999.;
472     Int_t iMultZNCTDC=0, iMultZPCTDC=0, iMultZEM1TDC=0, iMultZEM2TDC=0, iMultZNATDC=0, iMultZPATDC=0;
473     
474     const Float_t x[4] = {-1.75, 1.75, -1.75, 1.75};
475     const Float_t y[4] = {-1.75, -1.75, 1.75, 1.75};
476     const Float_t alpha=0.5;
477     Float_t numXZNC=0., numYZNC=0., denZNC=0., wZNC=0.; 
478     Float_t numXZNA=0., numYZNA=0., denZNA=0., wZNA=0.; 
479     
480     rawReader->Reset();
481     AliZDCRawStream stream(rawReader);
482     while(stream.Next()){
483
484       if(stream.IsADCDataWord() && 
485          (stream.GetADCModule()==0 || stream.GetADCModule()==1)){
486        
487          Int_t det = stream.GetSector(0);
488          Int_t quad = stream.GetSector(1);
489          Int_t gain = stream.GetADCGain();
490          Int_t pedindex=0;
491          
492          // Stuff for pedestal subtraction
493          if(quad != 5){ // ZDCs (not reference PTMs)
494           Float_t pedSubVal=-99.;
495           if(det == 1){    
496             pedindex = quad;
497             if(gain == 0){
498               pedSubVal = (Float_t) (stream.GetADCValue()-meanPed[pedindex]); 
499               zncSignal  += pedSubVal; 
500               isZNCFired = kTRUE;
501               if(quad!=0){
502                 zncSumQ += pedSubVal;
503                 if(pedSubVal>0.&& zncpmC>50.){
504                   wZNC = TMath::Power(pedSubVal, alpha);
505                   numXZNC += x[quad-1]*wZNC;
506                   numYZNC += y[quad-1]*wZNC;
507                   denZNC += wZNC;
508                 }
509               }
510               else{
511                 zncpmC = pedSubVal;
512                 GetRawsData(10)->Fill(zncpmC);
513               }
514               indZNC++;
515               
516               GetRawsData(20)->Fill(pedindex, pedSubVal);
517             }
518           }
519           else if(det == 2){ 
520             pedindex = quad+5;
521             if(gain == 0){
522               pedSubVal = (Float_t) (stream.GetADCValue()-meanPed[pedindex]); 
523               zpcSignal += pedSubVal; 
524               isZPCFired = kTRUE;
525               if(quad!=0) zpcSumQ += pedSubVal;
526               else{
527                 zpcpmC = pedSubVal;
528                 GetRawsData(12)->Fill(zpcpmC);
529               }
530               indZPC++;
531               
532               GetRawsData(20)->Fill(pedindex, pedSubVal);
533             }
534           }
535           else if(det == 3){ 
536             pedindex = quad+9;
537             if(quad==1){     
538               if(gain == 0){
539                 pedSubVal = (Float_t) (stream.GetADCValue()-meanPed[pedindex]); 
540                 GetRawsData(4)->Fill(pedSubVal);
541               
542                 GetRawsData(20)->Fill(pedindex, pedSubVal);
543               }
544             }
545             else if(quad==2){ 
546               if(gain == 0){
547                 pedSubVal = (Float_t) (stream.GetADCValue()-meanPed[pedindex]); 
548                 GetRawsData(5)->Fill(pedSubVal); 
549               
550                 GetRawsData(20)->Fill(pedindex, pedSubVal);
551               }
552             }
553           }
554           else if(det == 4){       
555             pedindex = quad+12;
556             if(gain == 0){
557               pedSubVal = (Float_t) (stream.GetADCValue()-meanPed[pedindex]); 
558               znaSignal  += pedSubVal; 
559               isZNAFired = kTRUE;
560               if(quad!=0){
561                 znaSumQ += pedSubVal;
562                 if(pedSubVal>0.&& znapmC>50.) {
563                   wZNA = TMath::Power(pedSubVal, alpha);
564                   numXZNA += x[quad-1]*wZNA;
565                   numYZNA += y[quad-1]*wZNA;
566                   denZNA += wZNA;
567                 }
568               }
569               else{
570                 znapmC = pedSubVal;
571                 GetRawsData(11)->Fill(znapmC);
572               }
573               indZNA++;
574               
575               GetRawsData(20)->Fill(pedindex, pedSubVal);
576               //GetRawsData(22)->SetBinContent(pedindex+1, GetRawsData(20)->GetBinContent(pedindex+1));
577               //GetRawsData(22)->SetBinError(pedindex+1, GetRawsData(20)->GetBinError(pedindex+1));
578             }
579           }
580           else if(det == 5){
581             pedindex = quad+17;
582             if(gain == 0){
583               pedSubVal = (Float_t) (stream.GetADCValue()-meanPed[pedindex]); 
584               zpaSignal  += pedSubVal; 
585               isZPAFired = kTRUE;
586               if(quad!=0) zpaSumQ += pedSubVal;
587               else{
588                 zpapmC = pedSubVal;
589                 GetRawsData(13)->Fill(zpapmC);
590               }
591               indZPA++;
592               
593               GetRawsData(20)->Fill(pedindex, pedSubVal);
594               //GetRawsData(22)->SetBinContent(pedindex+1, GetRawsData(20)->GetBinContent(pedindex+1));
595               //GetRawsData(22)->SetBinError(pedindex+1, GetRawsData(20)->GetBinError(pedindex+1));
596             }
597           }
598                          
599          }
600
601          if(isZNCFired && indZNC==5){
602            GetRawsData(0)->Fill(zncSignal);
603            GetRawsData(6)->Fill(zncpmC, zncSumQ);
604            GetRawsData(14)->Fill(zncSumQ); 
605            //
606            Float_t xZNC, yZNC;           
607            if(denZNC!=0){
608              xZNC = numXZNC/denZNC;
609              yZNC = numYZNC/denZNC;
610            } 
611            else xZNC = yZNC = 999.;
612            GetRawsData(24)->Fill(xZNC, yZNC);
613          }
614          if(isZPCFired && indZPC==5){
615            GetRawsData(2)->Fill(zpcSignal);
616            GetRawsData(8)->Fill(zpcpmC, zpcSumQ);
617            GetRawsData(16)->Fill(zpcSumQ); 
618          }
619          if(isZNAFired && indZNA==5){ 
620            GetRawsData(1)->Fill(znaSignal);
621            GetRawsData(7)->Fill(znapmC, znaSumQ);
622            GetRawsData(15)->Fill(znaSumQ); 
623            //
624            Float_t xZNA, yZNA;
625            if(denZNA!=0){
626              xZNA = numXZNA/denZNA;
627              yZNA = numYZNA/denZNA;
628            } 
629            else xZNA = yZNA = 999.;
630            GetRawsData(25)->Fill(xZNA, yZNA);
631          }
632          if(isZPAFired && indZPA==5){ 
633            GetRawsData(3)->Fill(zpaSignal);
634            GetRawsData(9)->Fill(zpapmC, zpaSumQ);
635            GetRawsData(17)->Fill(zpaSumQ); 
636          }
637          
638          if(indZNC==5){
639            zncSignal = zncSumQ = zncpmC = 0;
640            isZNCFired=kFALSE; indZNC=0;
641          }
642          if(indZPC==5){
643            zpcSignal = zpcSumQ = zpcpmC = 0;
644            isZPCFired=kFALSE; indZPC=0;
645          }
646          if(indZNA==5){
647            znaSignal = znaSumQ = znapmC = 0;
648            isZNAFired=kFALSE; indZNA=0;
649          }
650          if(indZPA==5){
651            zpaSignal = zpaSumQ = zpapmC = 0;
652            isZPAFired=kFALSE; indZPA=0;
653          }
654          
655       } //IsADCDataWord && signal ADCs
656       else if(stream.IsZDCTDCDatum()){
657          if(stream.GetChannel()==1){
658             zncTDC[iMultZNCTDC] = (0.025*(stream.GetZDCTDCDatum()));
659             iMultZNCTDC++;
660          }
661          else if(stream.GetChannel()==3){
662             zpcTDC[iMultZPCTDC] = (0.025*(stream.GetZDCTDCDatum()));
663             iMultZPCTDC++;
664          }
665          else if(stream.GetChannel()==5){
666             znaTDC[iMultZNATDC] = (0.025*(stream.GetZDCTDCDatum()));
667             iMultZNATDC++;
668          }
669          else if(stream.GetChannel()==7){
670             zpaTDC[iMultZPATDC] = (0.025*(stream.GetZDCTDCDatum()));
671             iMultZPATDC++;
672          }
673          else if(stream.GetChannel()==8){
674             zem1TDC[iMultZEM1TDC] = (0.025*(stream.GetZDCTDCDatum()));
675             iMultZEM1TDC++;
676          }
677          else if(stream.GetChannel()==9){
678             zem2TDC[iMultZEM2TDC] = (0.025*(stream.GetZDCTDCDatum()));
679             iMultZEM2TDC++;
680          }
681          else if(stream.GetChannel()==15) tdcGate = (0.025*(stream.GetZDCTDCDatum()));
682          
683          if(stream.GetChannel()==16 && tdcGate!=-999.){
684            for(Int_t iHit=0; iHit<10; iHit++){
685               if(zncTDC[iHit]!=-999.){
686                 GetRawsData(21)->Fill(1, zncTDC[iHit]-tdcGate);
687                 //GetRawsData(23)->SetBinContent(1, GetRawsData(21)->GetBinContent(1));
688               }
689               if(zpcTDC[iHit]!=-999.){
690                 Float_t diffZPC = zpcTDC[iHit]-tdcGate;
691                 GetRawsData(19)->Fill(diffZPC);
692                 GetRawsData(21)->Fill(2, diffZPC);
693                 //GetRawsData(23)->SetBinContent(2,  GetRawsData(21)->GetBinContent(3));
694               }
695               if(znaTDC[iHit]!=-999.){
696                 GetRawsData(21)->Fill(3, znaTDC[iHit]-tdcGate);
697                 //GetRawsData(23)->SetBinContent(3,  GetRawsData(21)->GetBinContent(5));
698               }
699               if(zpaTDC[iHit]!=-999.){
700                 GetRawsData(21)->Fill(4, zpaTDC[iHit]-tdcGate);
701                 //GetRawsData(23)->SetBinContent(4,  GetRawsData(21)->GetBinContent(7));
702               }
703               if(zem1TDC[iHit]!=-999.){
704                 Float_t diffZEM1 = zem1TDC[iHit]-tdcGate;
705                 GetRawsData(18)->Fill(diffZEM1);
706                 GetRawsData(21)->Fill(5, diffZEM1);
707                 //GetRawsData(23)->SetBinContent(5,  GetRawsData(21)->GetBinContent(8));
708               }
709               if(zem2TDC[iHit]!=-999.){
710                 GetRawsData(21)->Fill(6, zem2TDC[iHit]-tdcGate);
711                 //GetRawsData(23)->SetBinContent(6,  GetRawsData(21)->GetBinContent(9));
712               }
713            }
714            //
715            tdcGate = -999.;
716            for(Int_t i=0; i<10; i++){
717               zpcTDC[i] = zem1TDC[i] = -999.;
718            } 
719          }
720       }
721     
722     } //stream.Next()
723 //  } // check on event type
724 //  else{
725 //    AliDebug(1,Form("Skipping non-physics event for QA -> event type %d \n", rawReader->GetType())); 
726 //  }
727 }
728
729 //____________________________________________________________________________
730 void AliZDCQADataMakerRec::MakeRecPoints(TTree * clustersTree)
731 {
732   // Filling QA histos from RecPoints
733
734   TBranch *branch = clustersTree->GetBranch("ZDC");
735   if(!branch){ 
736     AliError("Can't get the ZDC branch for rec points!");
737     return;
738   }
739   
740   if(!GetRecPointsData(0)) InitRecPoints() ;
741
742   Float_t sum_ZNC=0., sum_ZNA=0., sum_ZPC=0., sum_ZPA=0.;
743   Float_t sumQ_ZNC=0., sumQ_ZNA=0., sumQ_ZPC=0., sumQ_ZPA=0.;
744
745   AliZDCReco reco;
746   AliZDCReco* preco = &reco;
747   clustersTree->SetBranchAddress("ZDC", &preco);
748
749   clustersTree->GetEntry(0);
750   for(Int_t i=0; i<5; i++){
751     sum_ZNC += reco.GetZN1HREnTow(i);
752     sum_ZPC += reco.GetZN2HREnTow(i);
753     sum_ZNA += reco.GetZP1HREnTow(i);
754     sum_ZPA += reco.GetZP2HREnTow(i);
755     if(i==0){
756       GetRecPointsData(8)->Fill(reco.GetZN1HREnTow(i));
757       GetRecPointsData(9)->Fill(reco.GetZN2HREnTow(i));
758       GetRecPointsData(10)->Fill(reco.GetZP1HREnTow(i));
759       GetRecPointsData(11)->Fill(reco.GetZP2HREnTow(i));
760     }
761     else{
762       sumQ_ZNC += reco.GetZN1HREnTow(i);
763       sumQ_ZPC += reco.GetZN2HREnTow(i);
764       sumQ_ZNA += reco.GetZP1HREnTow(i);
765       sumQ_ZPA += reco.GetZP2HREnTow(i);
766     }
767   }
768   
769   GetRecPointsData(0)->Fill(sum_ZNC);
770   GetRecPointsData(1)->Fill(sum_ZNA);
771   GetRecPointsData(2)->Fill(sum_ZPC);
772   GetRecPointsData(3)->Fill(sum_ZPA);
773   //
774   GetRecPointsData(4)->Fill(sumQ_ZNC);
775   GetRecPointsData(5)->Fill(sumQ_ZNA);
776   GetRecPointsData(6)->Fill(sumQ_ZPC);
777   GetRecPointsData(7)->Fill(sumQ_ZPA);
778   
779 }  
780
781 //____________________________________________________________________________
782 void AliZDCQADataMakerRec::MakeESDs(AliESDEvent * esd)
783 {
784   // make QA data from ESDs
785   //
786   
787   // Check id histograms already created for this Event Specie
788   if(!GetESDsData(0)) InitESDs() ;
789
790   AliESDZDC * zdcESD =  esd->GetESDZDC();
791   //
792   TString beamType = esd->GetBeamType();
793   Double_t centr_ZNC[2]={999.,999}, centr_ZNA[2]={999.,999};
794   if(((beamType.CompareTo("pp"))==0) || ((beamType.CompareTo("p-p"))==0)
795      ||((beamType.CompareTo("PP"))==0) || ((beamType.CompareTo("P-P"))==0)){
796     zdcESD->GetZNCentroidInpp(centr_ZNC, centr_ZNA);
797   }
798   else if((beamType.CompareTo("A-A")) == 0 || (beamType.CompareTo("Pb-Pb")) == 0){
799     Float_t beamEne = esd->GetBeamEnergy();
800     zdcESD->GetZNCentroidInPbPb(beamEne, centr_ZNC, centr_ZNA);
801   }
802   else printf("\n WARNING!!! AliZDCQADataMakerRec::MakeESDs: can't calculate centroids for beam type: %s\n\n",beamType.Data());
803   GetESDsData(0)->Fill(centr_ZNC[0], centr_ZNC[1]);
804   GetESDsData(1)->Fill(centr_ZNA[0], centr_ZNA[1]);
805
806   //
807   GetESDsData(2)->Fill(esd->GetZDCN1Energy());
808   GetESDsData(3)->Fill(esd->GetZDCN2Energy());
809   GetESDsData(4)->Fill(esd->GetZDCP1Energy());
810   GetESDsData(5)->Fill(esd->GetZDCP2Energy());
811   GetESDsData(6)->Fill(esd->GetZDCEMEnergy(0));
812   GetESDsData(7)->Fill(esd->GetZDCEMEnergy(1));
813   //
814   Double_t sumQZNC=0., sumQZPC=0., sumQZNA=0., sumQZPA=0.;
815   Double_t sumQZNC_lg=0., sumQZPC_lg=0., sumQZNA_lg=0., sumQZPA_lg=0.;
816   //
817   const Double_t *towZNC, *towZPC, *towZNA, *towZPA;
818   const Double_t *towZNC_lg, *towZPC_lg, *towZNA_lg, *towZPA_lg;
819   //
820   towZNC = zdcESD->GetZN1TowerEnergy();
821   towZPC = zdcESD->GetZP1TowerEnergy();
822   towZNA = zdcESD->GetZN2TowerEnergy();
823   towZPA = zdcESD->GetZP2TowerEnergy();
824   //
825   towZNC_lg = zdcESD->GetZN1TowerEnergyLR();
826   towZPC_lg = zdcESD->GetZP1TowerEnergyLR();
827   towZNA_lg = zdcESD->GetZN2TowerEnergyLR();
828   towZPA_lg = zdcESD->GetZP2TowerEnergyLR();
829   //
830   for(Int_t i=0; i<5; i++){
831      if(i==0){
832        GetESDsData(12)->Fill(towZNC[i]);
833        GetESDsData(13)->Fill(towZNA[i]);
834        GetESDsData(14)->Fill(towZPC[i]);
835        GetESDsData(15)->Fill(towZPA[i]);
836        //
837        GetESDsData(20)->Fill(towZNC_lg[i]);
838        GetESDsData(21)->Fill(towZNA_lg[i]);
839        GetESDsData(22)->Fill(towZPC_lg[i]);
840        GetESDsData(23)->Fill(towZPA_lg[i]);
841      }
842      else{
843        sumQZNC += towZNC[i];
844        sumQZPC += towZPC[i];
845        sumQZNA += towZNA[i];
846        sumQZPA += towZPA[i];
847        //
848        sumQZNC_lg += towZNC_lg[i];
849        sumQZPC_lg += towZPC_lg[i];
850        sumQZNA_lg += towZNA_lg[i];
851        sumQZPA_lg += towZPA_lg[i];
852      }
853   }
854   GetESDsData(8)->Fill(sumQZNC);
855   GetESDsData(9)->Fill(sumQZNA);
856   GetESDsData(10)->Fill(sumQZPC);
857   GetESDsData(11)->Fill(sumQZPA);
858   //
859   GetESDsData(16)->Fill(sumQZNC_lg);
860   GetESDsData(17)->Fill(sumQZNA_lg);
861   GetESDsData(18)->Fill(sumQZPC_lg);
862   GetESDsData(19)->Fill(sumQZPA_lg);
863 }
864
865 //____________________________________________________________________________
866 void AliZDCQADataMakerRec::StartOfDetectorCycle()
867 {
868   //Detector specific actions at start of cycle
869
870   fPedCalibData = GetPedCalibData();
871   
872 }
873
874 //____________________________________________________________________________ 
875 void AliZDCQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
876 {
877   //Detector specific actions at end of cycle
878   // do the QA checking
879   if( task == AliQAv1::kRAWS){
880      if (!GetRawsData(20) || !GetRawsData(21) || !GetRawsData(24)  || !GetRawsData(25)) {
881          printf("  WARNING!!! AliZDCQADataMaker Rec -> Not all histogram for DQM found!\n"); 
882      }
883      else{
884        for(Int_t ibin=1; ibin<=GetRawsData(20)->GetNbinsX(); ibin++){
885           GetRawsData(22)->SetBinContent(ibin, GetRawsData(20)->GetBinContent(ibin)); 
886           GetRawsData(22)->SetBinError(ibin, GetRawsData(20)->GetBinError(ibin));
887           GetRawsData(22)->SetLineColor(kBlue); GetRawsData(22)->SetLineWidth(2);
888        }
889        for(Int_t ibin=1; ibin<=GetRawsData(21)->GetNbinsX(); ibin++){
890           GetRawsData(23)->SetBinContent(ibin, GetRawsData(21)->GetBinContent(ibin)); 
891           GetRawsData(23)->SetBinError(ibin, GetRawsData(21)->GetBinError(ibin));
892           GetRawsData(23)->SetLineColor(kAzure-3); GetRawsData(23)->SetLineWidth(2);
893        }
894      }
895   }
896         
897   AliQAChecker::Instance()->Run(AliQAv1::kZDC, task, list) ;  
898 }