]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/AliZDCQADataMakerRec.cxx
Adding some features for DQM shifter histos
[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,0.,1200.);
178   TH1F * hZEM2Spectrum = new TH1F("hZEM2Spectrum","ZEM2 spectrum;Amplitude [ADC counts];Counts",100,0.,1200.);
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,7.,1407.,50,7., 1407.);
187   TH2F * hZPCpmCvsPMq = new TH2F("hZPCpmCvsPMq", "ZPC;PMC [ADC counts];Sum(PMQ) [ADC counts]",50,7.,1407.,50,7., 1407.);
188   TH2F * hZNApmCvsPMq = new TH2F("hZNApmCvsPMq", "ZNA;PMC [ADC counts];Sum(PMQ) [ADC counts]",50,7.,1407.,50,7., 1407.);
189   TH2F * hZPApmCvsPMq = new TH2F("hZPApmCvsPMq", "ZPA;PMC [ADC counts];Sum(PMQ) [ADC counts]",50,7.,1407.,50,7., 1407.);
190   Add2RawsList(hZNCpmCvsPMq, 6, !expert, image);
191   Add2RawsList(hZPCpmCvsPMq, 7, !expert, image);
192   Add2RawsList(hZNApmCvsPMq, 8, !expert, image);
193   Add2RawsList(hZPApmCvsPMq, 9, !expert, image);
194     
195   TH1F * hRawPMCZNC = new TH1F("hRawPMCZNC", "Raw ZNC PMC;Amplitude [ADC counts];Counts",100, 0., 1200.);
196   TH1F * hRawPMCZNA = new TH1F("hRawPMCZNA", "Raw ZNA PMC;Amplitude [ADC counts];Counts",100, 0., 1200.);
197   TH1F * hRawPMCZPC = new TH1F("hRawPMCZPC", "Raw ZPC PMC;Amplitude [ADC counts];Counts",100, 0., 1200.);
198   TH1F * hRawPMCZPA = new TH1F("hRawPMCZPA", "Raw ZPA PMC;Amplitude [ADC counts];Counts",100, 0., 1200.);
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 }
213
214 //____________________________________________________________________________
215 void AliZDCQADataMakerRec::InitRecPoints()
216 {
217   // create Digits histograms in Digits subdir
218   const Bool_t expert = kTRUE ; 
219   const Bool_t image  = kTRUE ; 
220   //
221   // ------------------- HIGH GAIN CHAIN ---------------------------
222   TH1F * hRecZNCTot = new TH1F("hRecZNCTot", "Rec signal in ZNC;Amplitude [ADC counts];Counts", 100, 0., 2000.);
223   TH1F * hRecZNATot = new TH1F("hRecZNATot", "Rec signal in ZNA;Amplitude [ADC counts];Counts", 100, 0., 2000.);
224   TH1F * hRecZPCTot = new TH1F("hRecZPCTot", "Rec signal in ZPC;Amplitude [ADC counts];Counts", 100, 0., 8000.);
225   TH1F * hRecZPATot = new TH1F("hRecZPATot", "Rec signal in ZPA;Amplitude [ADC counts];Counts", 100, 0., 8000.);
226   Add2RecPointsList(hRecZNCTot, 0, !expert, image);
227   Add2RecPointsList(hRecZNATot, 1, !expert, image);
228   Add2RecPointsList(hRecZPCTot, 2, !expert, image);
229   Add2RecPointsList(hRecZPATot, 3, !expert, image);
230   //
231   TH1F * hRecSumQZNC = new TH1F("hRecSumQZNC", "Rec summed 4 ZNC quadrants;Amplitude [ADC counts];Counts",100, 0., 2000.);
232   TH1F * hRecSumQZNA = new TH1F("hRecSumQZNA", "Rec summed 4 ZNA quadrants;Amplitude [ADC counts];Counts",100, 0., 2000.);
233   TH1F * hRecSumQZPC = new TH1F("hRecSumQZPC", "Rec summed 4 ZPC quadrants;Amplitude [ADC counts];Counts",100, 0., 2000.);
234   TH1F * hRecSumQZPA = new TH1F("hRecSumQZPA", "Rec summed 4 ZPA quadrants;Amplitude [ADC counts];Counts",100, 0., 2000.);
235   Add2RecPointsList(hRecSumQZNC, 4, expert, !image);
236   Add2RecPointsList(hRecSumQZNA, 5, expert, !image);
237   Add2RecPointsList(hRecSumQZPC, 6, expert, !image);
238   Add2RecPointsList(hRecSumQZPA, 7, expert, !image);
239   //
240   TH1F * hRecPMCZNC = new TH1F("hRecPMCZNC", "Rec common ZNC PMT;Amplitude [ADC counts];Counts",100, 0., 2000.);
241   TH1F * hRecPMCZNA = new TH1F("hRecPMCZNA", "Rec common ZNA PMT;Amplitude [ADC counts];Counts",100, 0., 2000.);
242   TH1F * hRecPMCZPC = new TH1F("hRecPMCZPC", "Rec common ZPC PMT;Amplitude [ADC counts];Counts",100, 0., 2000.);
243   TH1F * hRecPMCZPA = new TH1F("hRecPMCZPA", "Rec common ZPA PMT;Amplitude [ADC counts];Counts",100, 0., 2000.);
244   Add2RecPointsList(hRecPMCZNC, 8 , expert, !image);
245   Add2RecPointsList(hRecPMCZNA, 9 , expert, !image);
246   Add2RecPointsList(hRecPMCZPC, 10, expert, !image);
247   Add2RecPointsList(hRecPMCZPA, 11, expert, !image); 
248 }
249
250 //____________________________________________________________________________
251 void AliZDCQADataMakerRec::InitESDs()
252 {
253   //Booking ESDs histograms
254   //
255   const Bool_t expert = kTRUE ; 
256   const Bool_t image  = kTRUE ; 
257   
258   TH2F * hZNCcentr  = new TH2F("hZNCcentr", "Centroid in ZNC;X [cm];Y[cm]", 100, -5.,5.,100,-5.,5.);
259   TH2F * hZNAcentr  = new TH2F("hZNAcentr", "Centroid in ZNA;X [cm];Y[cm]", 100, -5.,5.,100,-5.,5.);
260   Add2ESDsList(hZNCcentr, 0, !expert, image);
261   Add2ESDsList(hZNAcentr, 1, !expert, image);
262   //
263   // ------------------- HIGH GAIN CHAIN ---------------------------
264   TH1F * hESDZNCTot = new TH1F("hESDZNCTot", "Energy in ZNC", 100, 0., 4000.);
265   TH1F * hESDZNATot = new TH1F("hESDZNATot", "Energy in ZNA", 100, 0., 4000.);
266   TH1F * hESDZPCTot = new TH1F("hESDZPCTot", "Energy in ZPC", 100, 0., 4000.);
267   TH1F * hESDZPATot = new TH1F("hESDZPATot", "Energy in ZPA", 100, 0., 4000.);
268   Add2ESDsList(hESDZNCTot, 2, !expert, image);
269   Add2ESDsList(hESDZNATot, 3, !expert, image);
270   Add2ESDsList(hESDZPCTot, 4, !expert, image);
271   Add2ESDsList(hESDZPATot, 5, !expert, image);
272   //
273   TH1F * hESDZEM1 = new TH1F("hESDZEM1", "Energy in ZEM1", 100, 0., 2000.);
274   TH1F * hESDZEM2 = new TH1F("hESDZEM2", "Energy in ZEM2", 100, 0., 2000.);
275   Add2ESDsList(hESDZEM1,6, !expert, image);
276   Add2ESDsList(hESDZEM2,7, !expert, image);
277   //
278   TH1F * hESDSumQZNC = new TH1F("hESDSumQZNC", "Sum of 4 ZNC energy",100, 0., 2000.);
279   TH1F * hESDSumQZNA = new TH1F("hESDSumQZNA", "Sum of 4 ZNA energy",100, 0., 2000.);
280   TH1F * hESDSumQZPC = new TH1F("hESDSumQZPC", "Sum of 4 ZPC energy",100, 0., 2000.);
281   TH1F * hESDSumQZPA = new TH1F("hESDSumQZPA", "Sum of 4 ZPA energy",100, 0., 2000.);
282   Add2ESDsList(hESDSumQZNC, 8, expert, !image);
283   Add2ESDsList(hESDSumQZNA, 9, expert, !image);
284   Add2ESDsList(hESDSumQZPC, 10, expert, !image);
285   Add2ESDsList(hESDSumQZPA, 11, expert, !image);
286   //
287   TH1F * hESDPMCZNC = new TH1F("hESDPMCZNC", "Energy in ZNC PMC",100, 0., 2000.);
288   TH1F * hESDPMCZNA = new TH1F("hESDPMCZNA", "Energy in ZNA PMC",100, 0., 2000.);
289   TH1F * hESDPMCZPC = new TH1F("hESDPMCZPC", "Energy in ZPC PMC",100, 0., 2000.);
290   TH1F * hESDPMCZPA = new TH1F("hESDPMCZPA", "Energy in ZPA PMC",100, 0., 2000.);
291   Add2ESDsList(hESDPMCZNC, 12, expert, !image);
292   Add2ESDsList(hESDPMCZNA, 13, expert, !image);
293   Add2ESDsList(hESDPMCZPC, 14, expert, !image);
294   Add2ESDsList(hESDPMCZPA, 15, expert, !image);
295   // 
296   // ------------------- LOW GAIN CHAIN ---------------------------
297   TH1F * hESDSumQZNClg = new TH1F("hESDSumQZNClg", "Sum of 4 lg ZNC sectors",100, 0., 4000.);
298   TH1F * hESDSumQZNAlg = new TH1F("hESDSumQZNAlg", "Sum of 4 lg ZNA sectors",100, 0., 4000.);
299   TH1F * hESDSumQZPClg = new TH1F("hESDSumQZPClg", "Sum of 4 lg ZPC sectors",100, 0., 4000.);
300   TH1F * hESDSumQZPAlg = new TH1F("hESDSumQZPAlg", "Sum of 4 lg ZPA sectors",100, 0., 4000.);
301   Add2ESDsList(hESDSumQZNClg, 16, expert, !image);
302   Add2ESDsList(hESDSumQZNAlg, 17, expert, !image);
303   Add2ESDsList(hESDSumQZPClg, 18, expert, !image);
304   Add2ESDsList(hESDSumQZPAlg, 19, expert, !image);
305   //
306   TH1F * hESDPMCZNClg = new TH1F("hESDPMCZNClg", "Signal in common ZNC lg PMT",100, 0., 4000.);
307   TH1F * hESDPMCZNAlg = new TH1F("hESDPMCZNAlg", "Signal in common ZNA lg PMT",100, 0., 4000.);
308   TH1F * hESDPMCZPClg = new TH1F("hESDPMCZPClg", "Signal in common ZPC lg PMT",100, 0., 4000.);
309   TH1F * hESDPMCZPAlg = new TH1F("hESDPMCZPAlg", "Signal in common ZPA lg PMT",100, 0., 4000.);
310   Add2ESDsList(hESDPMCZNClg, 20, expert, !image);
311   Add2ESDsList(hESDPMCZNAlg, 21, expert, !image);
312   Add2ESDsList(hESDPMCZPClg, 22, expert, !image);
313   Add2ESDsList(hESDPMCZPAlg, 23, expert, !image);
314 }
315
316 //___________________________________________________________________________
317 void AliZDCQADataMakerRec::MakeDigits(TTree *digitTree)
318 {
319   // makes data from Digit Tree
320   if(!GetDigitsData(0)) InitDigits();
321
322   if(!digitTree){
323     AliError("Can't get ZDC digit tree!!");
324     return; 
325   }     
326    
327   TBranch * branch = digitTree->GetBranch("ZDC");
328   if(!branch){
329     AliError("ZDC branch in digit tree not found"); 
330     return;
331   } 
332     
333   AliZDCDigit *digit = 0x0;
334   branch->SetAddress(&digit);
335      
336   Float_t adcSum_ZNC=0., adcSum_ZNA=0., adcSum_ZPC=0., adcSum_ZPA=0.;
337   Float_t adcSumQ_ZNC=0., adcSumQ_ZNA=0., adcSumQ_ZPC=0., adcSumQ_ZPA=0.;
338   Float_t adcSum_ZNC_lg=0., adcSum_ZNA_lg=0., adcSum_ZPC_lg=0., adcSum_ZPA_lg=0.;
339   Float_t adcSumQ_ZNC_lg=0., adcSumQ_ZNA_lg=0., adcSumQ_ZPC_lg=0., adcSumQ_ZPA_lg=0.;
340   
341   Int_t ndig = digitTree->GetEntries();
342   for(Int_t i=0; i<ndig; i++){
343       branch->GetEntry(i);
344       
345       if(digit->GetSector(0)==1 && digit->GetSector(1)!=5){
346           adcSum_ZNC += digit->GetADCValue(0);
347           adcSum_ZNC_lg += digit->GetADCValue(1);
348           //
349           if(digit->GetSector(1)!=0){
350               adcSumQ_ZNC += digit->GetADCValue(0);
351               adcSumQ_ZNC_lg+= digit->GetADCValue(1);
352           }
353           else{
354               GetDigitsData(8)->Fill(digit->GetADCValue(0));
355               GetDigitsData(20)->Fill(digit->GetADCValue(1));
356           }
357       }
358       else if(digit->GetSector(0)==2){
359           adcSum_ZPC += digit->GetADCValue(0);
360           adcSum_ZPC_lg += digit->GetADCValue(1);
361           //
362           if(digit->GetSector(1)!=0){
363               adcSumQ_ZPC += digit->GetADCValue(0);
364               adcSumQ_ZPC_lg+= digit->GetADCValue(1);
365           }
366           else{
367               GetDigitsData(10)->Fill(digit->GetADCValue(0));
368               GetDigitsData(22)->Fill(digit->GetADCValue(1));
369           }
370       }
371       else if(digit->GetSector(0)==4 && digit->GetSector(1)!=5){
372           adcSum_ZNA += digit->GetADCValue(0);
373           adcSum_ZNA_lg += digit->GetADCValue(1);
374           //
375           if(digit->GetSector(1)!=0){
376               adcSumQ_ZNA += digit->GetADCValue(0);
377               adcSumQ_ZNA_lg+= digit->GetADCValue(1);
378           }
379           else{
380               GetDigitsData(9)->Fill(digit->GetADCValue(0));
381               GetDigitsData(21)->Fill(digit->GetADCValue(1));
382           }
383       }
384       else if(digit->GetSector(0)==5){
385           adcSum_ZPA += digit->GetADCValue(0);
386           adcSum_ZPA_lg += digit->GetADCValue(1);
387           //
388           if(digit->GetSector(1)!=0){
389               adcSumQ_ZPA += digit->GetADCValue(0);
390               adcSumQ_ZPA_lg+= digit->GetADCValue(1);
391           }
392           else{
393               GetDigitsData(11)->Fill(digit->GetADCValue(0));
394               GetDigitsData(23)->Fill(digit->GetADCValue(1));
395           }
396       }
397   }
398   //
399   GetDigitsData(0)->Fill(adcSum_ZNC);
400   GetDigitsData(1)->Fill(adcSum_ZNA);
401   GetDigitsData(2)->Fill(adcSum_ZPC);
402   GetDigitsData(3)->Fill(adcSum_ZPA);
403   //
404   GetDigitsData(4)->Fill(adcSumQ_ZNC);
405   GetDigitsData(5)->Fill(adcSumQ_ZNA);
406   GetDigitsData(6)->Fill(adcSumQ_ZPC);
407   GetDigitsData(7)->Fill(adcSumQ_ZPA);
408   //
409   GetDigitsData(12)->Fill(adcSum_ZNC_lg);
410   GetDigitsData(13)->Fill(adcSum_ZNA_lg);
411   GetDigitsData(14)->Fill(adcSum_ZPC_lg);
412   GetDigitsData(15)->Fill(adcSum_ZPA_lg);
413   //
414   GetDigitsData(16)->Fill(adcSumQ_ZNC_lg);
415   GetDigitsData(17)->Fill(adcSumQ_ZNA_lg);
416   GetDigitsData(18)->Fill(adcSumQ_ZPC_lg);
417   GetDigitsData(19)->Fill(adcSumQ_ZPA_lg);
418   
419   delete digit;
420   digit=0;
421
422 }
423
424
425 //____________________________________________________________________________
426 void AliZDCQADataMakerRec::MakeRaws(AliRawReader *rawReader)
427 {
428   // Filling Raws QA histos
429   //
430   // Check if histograms already created for this Event Specie
431   if(!GetRawsData(0)) InitRaws();
432   
433   // Parameters for mean value pedestal subtraction
434   int const kNch = 24;
435   Float_t meanPed[2*kNch];    
436   for(Int_t jj=0; jj<2*kNch; jj++) meanPed[jj] = fPedCalibData->GetMeanPed(jj);
437   
438   AliZDCRawStream stream(rawReader);
439   while(stream.Next()){
440   
441     Float_t zncSignal=0., znaSignal=0., zpcSignal=0., zpaSignal=0.;
442     Float_t zncSumQ=0., znaSumQ=0., zpcSumQ=0., zpaSumQ=0.;
443     Float_t zncpmC=0., znapmC=0., zpcpmC=0., zpapmC=0.;
444     Bool_t isZNCFired=kFALSE, isZPCFired=kFALSE, isZNAFired=kFALSE, isZPAFired=kFALSE;
445     
446     if(stream.IsADCDataWord() && 
447      (stream.GetADCModule()==0 || stream.GetADCModule()==1)){
448      
449        Int_t det = stream.GetSector(0);
450        Int_t quad = stream.GetSector(1);
451        Int_t gain = stream.GetADCGain();
452        Int_t pedindex=0;
453        
454        // Stuff for pedestal subtraction
455        if(quad != 5){ // ZDCs (not reference PTMs)
456         if(det == 1){    
457           pedindex = quad;
458           if(gain == 0){
459             Float_t pedSubVal = (Float_t) (stream.GetADCValue()-meanPed[pedindex]); 
460             zncSignal  += pedSubVal; 
461             isZNCFired = kTRUE;
462             if(quad!=0) zncSumQ += pedSubVal;
463             else{
464               zncpmC = pedSubVal;
465               GetRawsData(10)->Fill(zncpmC);
466             }
467           }
468         }
469         else if(det == 2){ 
470           pedindex = quad+5;
471           if(gain == 0){
472             Float_t pedSubVal = (Float_t) (stream.GetADCValue()-meanPed[pedindex]); 
473             zpcSignal += pedSubVal; 
474             isZPCFired = kTRUE;
475             if(quad!=0) zpcSumQ += pedSubVal;
476             else{
477               zpcpmC = pedSubVal;
478               GetRawsData(12)->Fill(zpcpmC);
479             }
480           }
481         }
482         else if(det == 3){ 
483           pedindex = quad+9;
484           if(quad==1){     
485             if(gain == 0){
486               Float_t pedSubVal = (Float_t) (stream.GetADCValue()-meanPed[pedindex]); 
487               GetRawsData(4)->Fill(pedSubVal);
488             }
489           }
490           else if(quad==2){ 
491             if(gain == 0){
492               Float_t pedSubVal = (Float_t) (stream.GetADCValue()-meanPed[pedindex]); 
493               GetRawsData(5)->Fill(pedSubVal); 
494             }
495           }
496         }
497         else if(det == 4){       
498           pedindex = quad+12;
499           if(gain == 0){
500             Float_t pedSubVal = (Float_t) (stream.GetADCValue()-meanPed[pedindex]); 
501             znaSignal  += pedSubVal; 
502             isZNAFired = kTRUE;
503             if(quad!=0) znaSumQ += pedSubVal;
504             else{
505               znapmC = pedSubVal;
506               GetRawsData(11)->Fill(znapmC);
507             }
508           }
509         }
510         else if(det == 5){
511           pedindex = quad+17;
512           if(gain == 0){
513             Float_t pedSubVal = (Float_t) (stream.GetADCValue()-meanPed[pedindex]); 
514             zpaSignal  += pedSubVal; 
515             isZPAFired = kTRUE;
516             if(quad!=0) zpaSumQ += pedSubVal;
517             else{
518               zpapmC = pedSubVal;
519               GetRawsData(13)->Fill(zpapmC);
520             }
521           }
522         }
523        }
524        //
525       if(isZNCFired){
526         GetRawsData(0)->Fill(zncSignal);
527         GetRawsData(6)->Fill(zncSumQ, zncpmC);
528         GetRawsData(14)->Fill(zncSumQ); 
529       }
530       if(isZPCFired){
531         GetRawsData(2)->Fill(zpcSignal);
532         GetRawsData(8)->Fill(zpcSumQ, zpcpmC);
533         GetRawsData(15)->Fill(zpcSumQ); 
534       }
535       if(isZNAFired){ 
536         GetRawsData(1)->Fill(znaSignal);
537         GetRawsData(7)->Fill(znaSumQ, znapmC);
538         GetRawsData(16)->Fill(znaSumQ); 
539       }
540       if(isZPAFired){ 
541         GetRawsData(3)->Fill(zpaSignal);
542         GetRawsData(9)->Fill(znaSumQ, zpapmC);
543         GetRawsData(17)->Fill(znaSumQ); 
544       }
545        
546     } //IsADCDataWord && signal ADCs
547     
548   } //stream.Next()
549
550 //   stream.Delete();
551 }
552
553 //____________________________________________________________________________
554 void AliZDCQADataMakerRec::MakeRecPoints(TTree * clustersTree)
555 {
556   // Filling QA histos from RecPoints
557
558   TBranch *branch = clustersTree->GetBranch("ZDC");
559   if(!branch){ 
560     AliError("Can't get the ZDC branch for rec points!");
561     return;
562   }
563   
564   if(!GetRecPointsData(0)) InitRecPoints() ;
565
566   Float_t sum_ZNC=0., sum_ZNA=0., sum_ZPC=0., sum_ZPA=0.;
567   Float_t sumQ_ZNC=0., sumQ_ZNA=0., sumQ_ZPC=0., sumQ_ZPA=0.;
568
569   AliZDCReco reco;
570   AliZDCReco* preco = &reco;
571   clustersTree->SetBranchAddress("ZDC", &preco);
572
573   clustersTree->GetEntry(0);
574   for(Int_t i=0; i<5; i++){
575     sum_ZNC += reco.GetZN1HREnTow(i);
576     sum_ZPC += reco.GetZN2HREnTow(i);
577     sum_ZNA += reco.GetZP1HREnTow(i);
578     sum_ZPA += reco.GetZP2HREnTow(i);
579     if(i==0){
580       GetRecPointsData(8)->Fill(reco.GetZN1HREnTow(i));
581       GetRecPointsData(9)->Fill(reco.GetZN2HREnTow(i));
582       GetRecPointsData(10)->Fill(reco.GetZP1HREnTow(i));
583       GetRecPointsData(11)->Fill(reco.GetZP2HREnTow(i));
584     }
585     else{
586       sumQ_ZNC += reco.GetZN1HREnTow(i);
587       sumQ_ZPC += reco.GetZN2HREnTow(i);
588       sumQ_ZNA += reco.GetZP1HREnTow(i);
589       sumQ_ZPA += reco.GetZP2HREnTow(i);
590     }
591   }
592   
593   GetRecPointsData(0)->Fill(sum_ZNC);
594   GetRecPointsData(1)->Fill(sum_ZNA);
595   GetRecPointsData(2)->Fill(sum_ZPC);
596   GetRecPointsData(3)->Fill(sum_ZPA);
597   //
598   GetRecPointsData(4)->Fill(sumQ_ZNC);
599   GetRecPointsData(5)->Fill(sumQ_ZNA);
600   GetRecPointsData(6)->Fill(sumQ_ZPC);
601   GetRecPointsData(7)->Fill(sumQ_ZPA);
602   
603 }  
604
605 //____________________________________________________________________________
606 void AliZDCQADataMakerRec::MakeESDs(AliESDEvent * esd)
607 {
608   // make QA data from ESDs
609   //
610   
611   // Check id histograms already created for this Event Specie
612   if(!GetESDsData(0)) InitESDs() ;
613
614   AliESDZDC * zdcESD =  esd->GetESDZDC();
615   //
616   TString beamType = esd->GetBeamType();
617   Double_t centr_ZNC[2]={999.,999}, centr_ZNA[2]={999.,999};
618   if(((beamType.CompareTo("pp"))==0) || ((beamType.CompareTo("p-p"))==0)
619      ||((beamType.CompareTo("PP"))==0) || ((beamType.CompareTo("P-P"))==0)){
620     zdcESD->GetZNCentroidInpp(centr_ZNC, centr_ZNA);
621   }
622   else if((beamType.CompareTo("A-A")) == 0 || (beamType.CompareTo("Pb-Pb")) == 0){
623     Float_t beamEne = esd->GetBeamEnergy();
624     zdcESD->GetZNCentroidInPbPb(beamEne, centr_ZNC, centr_ZNA);
625   }
626   else printf(" AliZDCQADataMakerRec::MakeESDs: can't calculate centroids for beam type: %s\n\n",beamType.Data());
627   GetESDsData(0)->Fill(centr_ZNC[0], centr_ZNC[1]);
628   GetESDsData(1)->Fill(centr_ZNA[0], centr_ZNA[1]);
629
630   //
631   GetESDsData(2)->Fill(esd->GetZDCN1Energy());
632   GetESDsData(3)->Fill(esd->GetZDCN2Energy());
633   GetESDsData(4)->Fill(esd->GetZDCP1Energy());
634   GetESDsData(5)->Fill(esd->GetZDCP2Energy());
635   GetESDsData(6)->Fill(esd->GetZDCEMEnergy(0));
636   GetESDsData(7)->Fill(esd->GetZDCEMEnergy(1));
637   //
638   Double_t sumQZNC=0., sumQZPC=0., sumQZNA=0., sumQZPA=0.;
639   Double_t sumQZNC_lg=0., sumQZPC_lg=0., sumQZNA_lg=0., sumQZPA_lg=0.;
640   //
641   const Double_t *towZNC, *towZPC, *towZNA, *towZPA;
642   const Double_t *towZNC_lg, *towZPC_lg, *towZNA_lg, *towZPA_lg;
643   //
644   towZNC = zdcESD->GetZN1TowerEnergy();
645   towZPC = zdcESD->GetZP1TowerEnergy();
646   towZNA = zdcESD->GetZN2TowerEnergy();
647   towZPA = zdcESD->GetZP2TowerEnergy();
648   //
649   towZNC_lg = zdcESD->GetZN1TowerEnergyLR();
650   towZPC_lg = zdcESD->GetZP1TowerEnergyLR();
651   towZNA_lg = zdcESD->GetZN2TowerEnergyLR();
652   towZPA_lg = zdcESD->GetZP2TowerEnergyLR();
653   //
654   for(Int_t i=0; i<5; i++){
655      if(i==0){
656        GetESDsData(12)->Fill(towZNC[i]);
657        GetESDsData(13)->Fill(towZNA[i]);
658        GetESDsData(14)->Fill(towZPC[i]);
659        GetESDsData(15)->Fill(towZPA[i]);
660        //
661        GetESDsData(20)->Fill(towZNC_lg[i]);
662        GetESDsData(21)->Fill(towZNA_lg[i]);
663        GetESDsData(22)->Fill(towZPC_lg[i]);
664        GetESDsData(23)->Fill(towZPA_lg[i]);
665      }
666      else{
667        sumQZNC += towZNC[i];
668        sumQZPC += towZPC[i];
669        sumQZNA += towZNA[i];
670        sumQZPA += towZPA[i];
671        //
672        sumQZNC_lg += towZNC_lg[i];
673        sumQZPC_lg += towZPC_lg[i];
674        sumQZNA_lg += towZNA_lg[i];
675        sumQZPA_lg += towZPA_lg[i];
676      }
677   }
678   GetESDsData(8)->Fill(sumQZNC);
679   GetESDsData(9)->Fill(sumQZNA);
680   GetESDsData(10)->Fill(sumQZPC);
681   GetESDsData(11)->Fill(sumQZPA);
682   //
683   GetESDsData(16)->Fill(sumQZNC_lg);
684   GetESDsData(17)->Fill(sumQZNA_lg);
685   GetESDsData(18)->Fill(sumQZPC_lg);
686   GetESDsData(19)->Fill(sumQZPA_lg);
687 }
688
689 //____________________________________________________________________________
690 void AliZDCQADataMakerRec::StartOfDetectorCycle()
691 {
692   //Detector specific actions at start of cycle
693
694   fPedCalibData = GetPedCalibData();
695   
696 }
697
698 //____________________________________________________________________________ 
699 void AliZDCQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
700 {
701   //Detector specific actions at end of cycle
702   // do the QA checking
703   for(Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
704      if( !AliQAv1::Instance()->IsEventSpecieSet(specie) ) continue; 
705      
706     if (!GetRawsData(4) || !GetRawsData(5) || !GetRawsData(6) || !GetRawsData(7) || 
707         !GetRawsData(8) || !GetRawsData(9) || !GetRawsData(10) || !GetRawsData(11) || 
708         !GetRawsData(12) || !GetRawsData(13)) {
709          printf("  WARNING!!! AliZDCQADataMaker Rec -> No histogram for DQM found!\n") ; 
710          continue;
711      }
712      
713      TLine* diag = new TLine(7., 7., 1407., 1407.);
714      diag->SetLineColor(kRed);
715      diag->SetLineWidth(2);
716      
717      ((TH2F*)GetRawsData(6))->GetListOfFunctions()->Add(diag);
718      ((TH2F*)GetRawsData(7))->GetListOfFunctions()->Add(diag);
719      ((TH2F*)GetRawsData(8))->GetListOfFunctions()->Add(diag);
720      ((TH2F*)GetRawsData(9))->GetListOfFunctions()->Add(diag);
721
722      GetRawsData(6)->SetOption("colz");
723      GetRawsData(7)->SetOption("colz");
724      GetRawsData(8)->SetOption("colz");
725      GetRawsData(9)->SetOption("colz");  
726      
727      GetRawsData(4)->SetLineColor(kBlue+1);  GetRawsData(4)->SetLineWidth(2);
728      GetRawsData(5)->SetLineColor(kBlue+2);  GetRawsData(5)->SetLineWidth(2);
729      GetRawsData(10)->SetLineColor(kBlue+3); GetRawsData(10)->SetLineWidth(2);
730      GetRawsData(11)->SetLineColor(kBlue+4); GetRawsData(11)->SetLineWidth(2);
731      GetRawsData(12)->SetLineColor(kBlue+5); GetRawsData(12)->SetLineWidth(2);
732      GetRawsData(13)->SetLineColor(kBlue+6); GetRawsData(13)->SetLineWidth(2);
733      
734      GetRawsData(4)->SetDrawOption("LOGY");
735      GetRawsData(5)->SetDrawOption("LOGY");
736      GetRawsData(10)->SetDrawOption("LOGY");
737      GetRawsData(11)->SetDrawOption("LOGY");
738      GetRawsData(12)->SetDrawOption("LOGY");
739      GetRawsData(13)->SetDrawOption("LOGY");
740
741   }
742         
743   AliQAChecker::Instance()->Run(AliQAv1::kZDC, task, list) ;  
744 }
745