1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 ///////////////////////////////////////////////////////////////////////
18 // Produces the data needed to calculate the ZDC quality assurance. //
19 // QA objects are 1 & 2 Dimensional histograms. //
20 // author: C. Oppedisano //
22 ///////////////////////////////////////////////////////////////////////
24 // --- ROOT system ---
25 #include <TClonesArray.h>
30 #include <Riostream.h>
31 // --- Standard library ---
33 // --- AliRoot header files ---
35 #include "AliQAChecker.h"
36 #include "AliRawReader.h"
37 #include "AliZDCQADataMakerRec.h"
38 #include "AliZDCRawStream.h"
39 #include "AliZDCDigit.h"
40 #include "AliESDZDC.h"
41 #include "AliESDEvent.h"
43 ClassImp(AliZDCQADataMakerRec)
45 //____________________________________________________________________________
46 AliZDCQADataMakerRec::AliZDCQADataMakerRec() :
47 AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kZDC), "ZDC Quality Assurance Data Maker"),
53 //____________________________________________________________________________
54 AliZDCQADataMakerRec::AliZDCQADataMakerRec(const AliZDCQADataMakerRec& qadm) :
59 SetName((const char*)qadm.GetName());
60 SetTitle((const char*)qadm.GetTitle());
63 //__________________________________________________________________
64 AliZDCQADataMakerRec& AliZDCQADataMakerRec::operator = (const AliZDCQADataMakerRec& qadm )
67 this->~AliZDCQADataMakerRec();
68 new(this) AliZDCQADataMakerRec(qadm);
72 //____________________________________________________________________________
74 void AliZDCQADataMakerRec::InitRaws()
76 // create Digits histograms in Digits subdir
78 // ------------------- HIGH GAIN CHAIN ---------------------------
79 TH1F * hRawZNCTot = new TH1F("hRawZNCTot", "Raw signal in ZNC;Amplitude [ADC counts];Counts", 100, 0., 6000.);
80 TH1F * hRawZNATot = new TH1F("hRawZNATot", "Raw signal in ZNA;Amplitude [ADC counts];Counts", 100, 0., 6000.);
81 TH1F * hRawZPCTot = new TH1F("hRawZPCTot", "Raw signal in ZPC;Amplitude [ADC counts];Counts", 100, 0., 10000.);
82 TH1F * hRawZPATot = new TH1F("hRawZPATot", "Raw signal in ZPA;Amplitude [ADC counts];Counts", 100, 0., 10000.);
83 Add2RawsList(hRawZNCTot, 0);
84 Add2RawsList(hRawZNATot, 1);
85 Add2RawsList(hRawZPCTot, 2);
86 Add2RawsList(hRawZPATot, 3);
88 TH1F * hRawSumQZNC = new TH1F("hRawSumQZNC", "Raw summed 4 ZNC quadrants;Amplitude [ADC counts];Counts",100, 0., 4000.);
89 TH1F * hRawSumQZNA = new TH1F("hRawSumQZNA", "Raw summed 4 ZNA quadrants;Amplitude [ADC counts];Counts",100, 0., 4000.);
90 TH1F * hRawSumQZPC = new TH1F("hRawSumQZPC", "Raw summed 4 ZPC quadrants;Amplitude [ADC counts];Counts",100, 0., 4000.);
91 TH1F * hRawSumQZPA = new TH1F("hRawSumQZPA", "Raw summed 4 ZPA quadrants;Amplitude [ADC counts];Counts",100, 0., 4000.);
92 Add2RawsList(hRawSumQZNC, 4, kTRUE);
93 Add2RawsList(hRawSumQZNA, 5, kTRUE);
94 Add2RawsList(hRawSumQZPC, 6, kTRUE);
95 Add2RawsList(hRawSumQZPA, 7, kTRUE);
97 TH1F * hRawPMCZNC = new TH1F("hRawPMCZNC", "Raw common ZNC PMT;Amplitude [ADC counts];Counts",100, 0., 4000.);
98 TH1F * hRawPMCZNA = new TH1F("hRawPMCZNA", "Raw common ZNA PMT;Amplitude [ADC counts];Counts",100, 0., 4000.);
99 TH1F * hRawPMCZPC = new TH1F("hRawPMCZPC", "Raw common ZPC PMT;Amplitude [ADC counts];Counts",100, 0., 4000.);
100 TH1F * hRawPMCZPA = new TH1F("hRawPMCZPA", "Raw common ZPA PMT;Amplitude [ADC counts];Counts",100, 0., 4000.);
101 Add2RawsList(hRawPMCZNC, 8 , kTRUE);
102 Add2RawsList(hRawPMCZNA, 9 , kTRUE);
103 Add2RawsList(hRawPMCZPC, 10, kTRUE);
104 Add2RawsList(hRawPMCZPA, 11, kTRUE);
106 /* // ------------------- LOW GAIN CHAIN ---------------------------
107 TH1F * hRawZNCTotlg = new TH1F("hRawZNCTotlg", "Rawit lg signal in ZNC", 100, 0., 6000.);
108 TH1F * hRawZNATotlg = new TH1F("hRawZNATotlg", "Rawit lg signal in ZNA", 100, 0., 6000.);
109 TH1F * hRawZPCTotlg = new TH1F("hRawZPCTotlg", "Rawit lg signal in ZPC", 100, 0., 10000.);
110 TH1F * hRawZPATotlg = new TH1F("hRawZPATotlg", "Rawit lg signal in ZPA", 100, 0., 10000.);
111 Add2RawsList(hRawZNCTotlg, 12);
112 Add2RawsList(hRawZNATotlg, 13);
113 Add2RawsList(hRawZPCTotlg, 14);
114 Add2RawsList(hRawZPATotlg, 15);
116 TH1F * hRawSumQZNClg = new TH1F("hRawSumQZNClg", "Raw summed 4 lg ZNC quadrants",100, 0., 4000.);
117 TH1F * hRawSumQZNAlg = new TH1F("hRawSumQZNAlg", "Raw summed 4 lg ZNA quadrants",100, 0., 4000.);
118 TH1F * hRawSumQZPClg = new TH1F("hRawSumQZPClg", "Raw summed 4 lg ZPC quadrants",100, 0., 4000.);
119 TH1F * hRawSumQZPAlg = new TH1F("hRawSumQZPAlg", "Raw summed 4 lg ZPA quadrants",100, 0., 4000.);
120 Add2RawsList(hRawSumQZNClg, 16, kTRUE);
121 Add2RawsList(hRawSumQZNAlg, 17, kTRUE);
122 Add2RawsList(hRawSumQZPClg, 18, kTRUE);
123 Add2RawsList(hRawSumQZPAlg, 19, kTRUE);
125 TH1F * hRawPMCZNClg = new TH1F("hRawPMCZNClg", "Raw common lg ZNC PMT",100, 0., 4000.);
126 TH1F * hRawPMCZNAlg = new TH1F("hRawPMCZNAlg", "Raw common lg ZNA PMT",100, 0., 4000.);
127 TH1F * hRawPMCZPClg = new TH1F("hRawPMCZPClg", "Raw common lg ZPC PMT",100, 0., 4000.);
128 TH1F * hRawPMCZPAlg = new TH1F("hRawPMCZPAlg", "Raw common lg ZPA PMT",100, 0., 4000.);
129 Add2RawsList(hRawPMCZNClg, 20, kTRUE);
130 Add2RawsList(hRawPMCZNAlg, 21, kTRUE);
131 Add2RawsList(hRawPMCZPClg, 22, kTRUE);
132 Add2RawsList(hRawPMCZPAlg, 23, kTRUE);*/
135 //____________________________________________________________________________
136 void AliZDCQADataMakerRec::InitDigits()
138 // create Digits histograms in Digits subdir
140 const Bool_t expert = kTRUE ;
141 const Bool_t image = kTRUE ;
143 // ------------------- HIGH GAIN CHAIN ---------------------------
144 TH1F * hDigZNCTot = new TH1F("hDigZNCTot", "Signal in ZNC;Amplitude [ADC counts];Counts", 100, 0., 6000.);
145 TH1F * hDigZNATot = new TH1F("hDigZNATot", "Signal in ZNA;Amplitude [ADC counts];Counts", 100, 0., 6000.);
146 TH1F * hDigZPCTot = new TH1F("hDigZPCTot", "Signal in ZPC;Amplitude [ADC counts];Counts", 100, 0., 6000.);
147 TH1F * hDigZPATot = new TH1F("hDigZPATot", "Signal in ZPA;Amplitude [ADC counts];Counts", 100, 0., 6000.);
148 Add2DigitsList(hDigZNCTot, 0, !expert, image);
149 Add2DigitsList(hDigZNATot, 1, !expert, image);
150 Add2DigitsList(hDigZPCTot, 2, !expert, image);
151 Add2DigitsList(hDigZPATot, 3, !expert, image);
153 TH1F * hDigSumQZNC = new TH1F("hDigSumQZNC", "Signal in 4 ZNC PMQ;Amplitude [ADC counts];Counts",100, 0., 4000.);
154 TH1F * hDigSumQZNA = new TH1F("hDigSumQZNA", "Signal in 4 ZNA PMQ;Amplitude [ADC counts];Counts",100, 0., 4000.);
155 TH1F * hDigSumQZPC = new TH1F("hDigSumQZPC", "Signal in 4 ZPC PMQ;Amplitude [ADC counts];Counts",100, 0., 4000.);
156 TH1F * hDigSumQZPA = new TH1F("hDigSumQZPA", "Signal in 4 ZPA PMQ;Amplitude [ADC counts];Counts",100, 0., 4000.);
157 Add2DigitsList(hDigSumQZNC, 4, expert, !image);
158 Add2DigitsList(hDigSumQZNA, 5, expert, !image);
159 Add2DigitsList(hDigSumQZPC, 6, expert, !image);
160 Add2DigitsList(hDigSumQZPA, 7, expert, !image);
162 TH1F * hDigPMCZNC = new TH1F("hDigPMCZNC", "Signal in ZNC PMC;Amplitude [ADC counts];Counts",100, 0., 4000.);
163 TH1F * hDigPMCZNA = new TH1F("hDigPMCZNA", "Signal in ZNA PMC;Amplitude [ADC counts];Counts",100, 0., 4000.);
164 TH1F * hDigPMCZPC = new TH1F("hDigPMCZPC", "Signal in ZPC PMC;Amplitude [ADC counts];Counts",100, 0., 4000.);
165 TH1F * hDigPMCZPA = new TH1F("hDigPMCZPA", "Signal in ZPA PMC;Amplitude [ADC counts];Counts",100, 0., 4000.);
166 Add2DigitsList(hDigPMCZNC, 8, expert, !image);
167 Add2DigitsList(hDigPMCZNA, 9, expert, !image);
168 Add2DigitsList(hDigPMCZPC, 10, expert, !image);
169 Add2DigitsList(hDigPMCZPA, 11, expert, !image);
171 // ------------------- LOW GAIN CHAIN ---------------------------
172 /* TH1F * hDigZNCTotlg = new TH1F("hDigZNCTotlg", "Digit lg signal in ZNC", 100, 0., 6000.);
173 TH1F * hDigZNATotlg = new TH1F("hDigZNATotlg", "Digit lg signal in ZNA", 100, 0., 6000.);
174 TH1F * hDigZPCTotlg = new TH1F("hDigZPCTotlg", "Digit lg signal in ZPC", 100, 0., 6000.);
175 TH1F * hDigZPATotlg = new TH1F("hDigZPATotlg", "Digit lg signal in ZPA", 100, 0., 6000.);
176 Add2DigitsList(hDigZNCTotlg, 12, !expert, image);
177 Add2DigitsList(hDigZNATotlg, 13, !expert, image);
178 Add2DigitsList(hDigZPCTotlg, 14, !expert, image);
179 Add2DigitsList(hDigZPATotlg, 15, !expert, image);
181 TH1F * hDigSumQZNClg = new TH1F("hDigSumQZNClg", "Signal in 4 ZNC PMQlg",100, 0., 4000.);
182 TH1F * hDigSumQZNAlg = new TH1F("hDigSumQZNAlg", "Signal in 4 ZNA PMQlg",100, 0., 4000.);
183 TH1F * hDigSumQZPClg = new TH1F("hDigSumQZPClg", "Signal in 4 ZPC PMQlg",100, 0., 4000.);
184 TH1F * hDigSumQZPAlg = new TH1F("hDigSumQZPAlg", "Signal in 4 ZPA PMQlg",100, 0., 4000.);
185 Add2DigitsList(hDigSumQZNClg, 16, expert, !image);
186 Add2DigitsList(hDigSumQZNAlg, 17, expert, !image);
187 Add2DigitsList(hDigSumQZPClg, 18, expert, !image);
188 Add2DigitsList(hDigSumQZPAlg, 19, expert, !image);
190 TH1F * hDigPMCZNClg = new TH1F("hDigPMCZNClg", "Signal in ZNC PMClg",100, 0., 4000.);
191 TH1F * hDigPMCZNAlg = new TH1F("hDigPMCZNAlg", "Signal in ZNA PMClg",100, 0., 4000.);
192 TH1F * hDigPMCZPClg = new TH1F("hDigPMCZPClg", "Signal in ZPC PMClg",100, 0., 4000.);
193 TH1F * hDigPMCZPAlg = new TH1F("hDigPMCZPAlg", "Signal in ZPA PMClg",100, 0., 4000.);
194 Add2DigitsList(hDigPMCZNClg, 20, expert, !image);
195 Add2DigitsList(hDigPMCZNAlg, 21, expert, !image);
196 Add2DigitsList(hDigPMCZPClg, 22, expert, !image);
197 Add2DigitsList(hDigPMCZPAlg, 23, expert, !image);
201 //____________________________________________________________________________
202 void AliZDCQADataMakerRec::InitESDs()
204 //Booking ESDs histograms
206 const Bool_t expert = kTRUE ;
207 const Bool_t image = kTRUE ;
209 TH2F * hZNC = new TH2F("hZNC", "Centroid in ZNC", 100, -5.,5.,100,-5.,5.);
210 TH2F * hZNA = new TH2F("hZNA", "Centroid in ZNA", 100, -5.,5.,100,-5.,5.);
211 Add2ESDsList(hZNC, 0, !expert, image);
212 Add2ESDsList(hZNA, 1, !expert, image);
214 // ------------------- HIGH GAIN CHAIN ---------------------------
215 TH1F * hESDZNCTot = new TH1F("hESDZNCTot", "Energy in ZNC", 100, 0., 6000.);
216 TH1F * hESDZNATot = new TH1F("hESDZNATot", "Energy in ZNA", 100, 0., 6000.);
217 TH1F * hESDZPCTot = new TH1F("hESDZPCTot", "Energy in ZPC", 100, 0., 10000.);
218 TH1F * hESDZPATot = new TH1F("hESDZPATot", "Energy in ZPA", 100, 0., 10000.);
219 Add2ESDsList(hESDZNCTot, 2, !expert, image);
220 Add2ESDsList(hESDZNATot, 3, !expert, image);
221 Add2ESDsList(hESDZPCTot, 4, !expert, image);
222 Add2ESDsList(hESDZPATot, 5, !expert, image);
224 TH1F * hESDSumQZNC = new TH1F("hESDSumQZNC", "Sum of 4 ZNC energy",100, 0., 4000.);
225 TH1F * hESDSumQZNA = new TH1F("hESDSumQZNA", "Sum of 4 ZNA energy",100, 0., 4000.);
226 TH1F * hESDSumQZPC = new TH1F("hESDSumQZPC", "Sum of 4 ZPC energy",100, 0., 4000.);
227 TH1F * hESDSumQZPA = new TH1F("hESDSumQZPA", "Sum of 4 ZPA energy",100, 0., 4000.);
228 Add2ESDsList(hESDSumQZNC, 6, expert, !image);
229 Add2ESDsList(hESDSumQZNA, 7, expert, !image);
230 Add2ESDsList(hESDSumQZPC, 8, expert, !image);
231 Add2ESDsList(hESDSumQZPA, 9, expert, !image);
233 TH1F * hESDPMCZNC = new TH1F("hESDPMCZNC", "Energy in common ZNC PMT",100, 0., 4000.);
234 TH1F * hESDPMCZNA = new TH1F("hESDPMCZNA", "Energy in common ZNA PMT",100, 0., 4000.);
235 TH1F * hESDPMCZPC = new TH1F("hESDPMCZPC", "Energy in common ZPC PMT",100, 0., 4000.);
236 TH1F * hESDPMCZPA = new TH1F("hESDPMCZPA", "Energy in common ZPA PMT",100, 0., 4000.);
237 Add2ESDsList(hESDPMCZNC, 10, expert, !image);
238 Add2ESDsList(hESDPMCZNA, 11, expert, !image);
239 Add2ESDsList(hESDPMCZPC, 12, expert, !image);
240 Add2ESDsList(hESDPMCZPA, 13, expert, !image);
242 /* // ------------------- LOW GAIN CHAIN ---------------------------
243 TH1F * hESDZNCTotlg = new TH1F("hESDZNCTotlg", "ESD lg signal in ZNC", 100, 0., 6000.);
244 TH1F * hESDZNATotlg = new TH1F("hESDZNATotlg", "ESD lg signal in ZNA", 100, 0., 6000.);
245 TH1F * hESDZPCTotlg = new TH1F("hESDZPCTotlg", "ESD lg signal in ZPC", 100, 0., 10000.);
246 TH1F * hESDZPATotlg = new TH1F("hESDZPATotlg", "ESD lg signal in ZPA", 100, 0., 10000.);
247 Add2ESDsList(hESDZNCTotlg, !expert, image);
248 Add2ESDsList(hESDZNATotlg, !expert, image);
249 Add2ESDsList(hESDZPCTotlg, !expert, image);
250 Add2ESDsList(hESDZPATotlg, !expert, image);
252 TH1F * hESDSumQZNClg = new TH1F("hESDSumQZNClg", "Sum of 4 lg ZNC sectors",100, 0., 4000.);
253 TH1F * hESDSumQZNAlg = new TH1F("hESDSumQZNAlg", "Sum of 4 lg ZNA sectors",100, 0., 4000.);
254 TH1F * hESDSumQZPClg = new TH1F("hESDSumQZPClg", "Sum of 4 lg ZPC sectors",100, 0., 4000.);
255 TH1F * hESDSumQZPAlg = new TH1F("hESDSumQZPAlg", "Sum of 4 lg ZPA sectors",100, 0., 4000.);
256 Add2ESDsList(hESDSumQZNClg, 18, expert, !image);
257 Add2ESDsList(hESDSumQZNAlg, 19, expert, !image);
258 Add2ESDsList(hESDSumQZPClg, 20, expert, !image);
259 Add2ESDsList(hESDSumQZPAlg, 21, expert, !image);
261 TH1F * hESDPMCZNClg = new TH1F("hESDPMCZNClg", "Signal in common ZNC lg PMT",100, 0., 4000.);
262 TH1F * hESDPMCZNAlg = new TH1F("hESDPMCZNAlg", "Signal in common ZNA lg PMT",100, 0., 4000.);
263 TH1F * hESDPMCZPClg = new TH1F("hESDPMCZPClg", "Signal in common ZPC lg PMT",100, 0., 4000.);
264 TH1F * hESDPMCZPAlg = new TH1F("hESDPMCZPAlg", "Signal in common ZPA lg PMT",100, 0., 4000.);
265 Add2ESDsList(hESDPMCZNClg, 22, expert, !image);
266 Add2ESDsList(hESDPMCZNAlg, 23, expert, !image);
267 Add2ESDsList(hESDPMCZPClg, 24, expert, !image);
268 Add2ESDsList(hESDPMCZPAlg, 25, expert, !image);*/
271 //____________________________________________________________________________
273 void AliZDCQADataMakerRec::MakeRaws(AliRawReader *rawReader)
275 // Filling Raws QA histos
277 // Check id histograms already created for this Event Specie
278 if ( ! GetRawsData(0) )
282 Float_t sum_ZNC=0., sum_ZNA=0., sum_ZPC=0., sum_ZPA=0.;
283 Float_t sumQ_ZNC=0., sumQ_ZNA=0., sumQ_ZPC=0., sumQ_ZPA=0.;
284 //Float_t sum_ZNC_lg=0., sum_ZNA_lg=0., sum_ZPC_lg=0., sum_ZPA_lg=0.;
285 //Float_t sumQ_ZNC_lg=0., sumQ_ZNA_lg=0., sumQ_ZPC_lg=0., sumQ_ZPA_lg=0.;
287 AliZDCRawStream stream(rawReader);
288 while(stream.Next()){
289 if(stream.IsADCDataWord() &&
290 (stream.GetADCModule()==0 || stream.GetADCModule()==1)){
291 if(stream.GetSector(0)==1){
292 if(stream.GetADCGain()==0){
293 sum_ZNC += stream.GetADCValue();
294 if(stream.GetSector(1)!=0) sumQ_ZNC += stream.GetADCValue();
295 else GetRawsData(8)->Fill(stream.GetADCValue());
298 sum_ZNC_lg += stream.GetADCValue();
299 if(stream.GetSector(1)!=0) sumQ_ZNC_lg += stream.GetADCValue();
300 else GetRawsData(20)->Fill(stream.GetADCValue());
303 else if(stream.GetSector(0)==2){
304 if(stream.GetADCGain()==0){
305 sum_ZPC += stream.GetADCValue();
306 if(stream.GetSector(1)!=0) sumQ_ZPC += stream.GetADCValue();
307 else GetRawsData(10)->Fill(stream.GetADCValue());
310 sum_ZPC_lg += stream.GetADCValue();
311 if(stream.GetSector(1)!=0) sumQ_ZPC_lg += stream.GetADCValue();
312 else GetRawsData(22)->Fill(stream.GetADCValue());
315 else if(stream.GetSector(0)==4){
316 if(stream.GetADCGain()==0){
317 sum_ZNA += stream.GetADCValue();
318 if(stream.GetSector(1)!=0) sumQ_ZNA += stream.GetADCValue();
319 else GetRawsData(9)->Fill(stream.GetADCValue());
322 sum_ZNA_lg += stream.GetADCValue();
323 if(stream.GetSector(1)!=0) sumQ_ZNA_lg += stream.GetADCValue();
324 else GetRawsData(21)->Fill(stream.GetADCValue());
327 else if(stream.GetSector(0)==5){
328 if(stream.GetADCGain()==0){
329 sum_ZPA += stream.GetADCValue();
330 if(stream.GetSector(1)!=0) sumQ_ZPA += stream.GetADCValue();
331 else GetRawsData(11)->Fill(stream.GetADCValue());
334 sum_ZPA_lg += stream.GetADCValue();
335 if(stream.GetSector(1)!=0) sumQ_ZPA_lg += stream.GetADCValue();
336 else GetRawsData(23)->Fill(stream.GetADCValue());
342 GetRawsData(0)->Fill(sum_ZNC);
343 GetRawsData(1)->Fill(sum_ZNA);
344 GetRawsData(2)->Fill(sum_ZPC);
345 GetRawsData(3)->Fill(sum_ZPA);
347 GetRawsData(4)->Fill(sumQ_ZNC);
348 GetRawsData(5)->Fill(sumQ_ZNA);
349 GetRawsData(6)->Fill(sumQ_ZPC);
350 GetRawsData(7)->Fill(sumQ_ZPA);
352 /*GetRawsData(12)->Fill(sum_ZNC_lg);
353 GetRawsData(13)->Fill(sum_ZNA_lg);
354 GetRawsData(14)->Fill(sum_ZPC_lg);
355 GetRawsData(15)->Fill(sum_ZPA_lg);
357 GetRawsData(16)->Fill(sumQ_ZNC_lg);
358 GetRawsData(17)->Fill(sumQ_ZNA_lg);
359 GetRawsData(18)->Fill(sumQ_ZPC_lg);
360 GetRawsData(19)->Fill(sumQ_ZPA_lg);*/
365 //___________________________________________________________________________
366 void AliZDCQADataMakerRec::MakeDigits(TTree *digitTree )
368 TBranch * branch = digitTree->GetBranch("ZDC");
370 AliError("ZDC branch in Digit Tree not found");
374 // Check id histograms already created for this Event Specie
375 if ( ! GetDigitsData(0) )
378 branch->SetAddress(&fDigit);
380 Int_t ndig = digitTree->GetEntries();
382 Float_t adcSum_ZNC=0., adcSum_ZNA=0., adcSum_ZPC=0., adcSum_ZPA=0.;
383 Float_t adcSumQ_ZNC=0., adcSumQ_ZNA=0., adcSumQ_ZPC=0., adcSumQ_ZPA=0.;
384 //Float_t adcSum_ZNC_lg=0., adcSum_ZNA_lg=0., adcSum_ZPC_lg=0., adcSum_ZPA_lg=0.;
385 //Float_t adcSumQ_ZNC_lg=0., adcSumQ_ZNA_lg=0., adcSumQ_ZPC_lg=0., adcSumQ_ZPA_lg=0.;
387 for(Int_t i = 0; i < ndig; i++){
388 digitTree->GetEntry(i);
389 if(fDigit->GetSector(0)==1){
390 adcSum_ZNC += fDigit->GetADCValue(0);
391 //adcSum_ZNC_lg += fDigit->GetADCValue(1);
393 if(fDigit->GetSector(1)!=0){
394 adcSumQ_ZNC += fDigit->GetADCValue(0);
395 //adcSumQ_ZNC_lg+= fDigit->GetADCValue(1);
398 GetDigitsData(8)->Fill(fDigit->GetADCValue(0));
399 //GetDigitsData(20)->Fill(fDigit->GetADCValue(1));
402 else if(fDigit->GetSector(0)==2){
403 adcSum_ZPC += fDigit->GetADCValue(0);
404 //adcSum_ZPC_lg += fDigit->GetADCValue(1);
406 if(fDigit->GetSector(1)!=0){
407 adcSumQ_ZPC += fDigit->GetADCValue(0);
408 //adcSumQ_ZPC_lg+= fDigit->GetADCValue(1);
411 GetDigitsData(10)->Fill(fDigit->GetADCValue(0));
412 //GetDigitsData(22)->Fill(fDigit->GetADCValue(1));
415 else if(fDigit->GetSector(0)==4){
416 adcSum_ZNA += fDigit->GetADCValue(0);
417 //adcSum_ZNA_lg += fDigit->GetADCValue(1);
419 if(fDigit->GetSector(1)!=0){
420 adcSumQ_ZNA += fDigit->GetADCValue(0);
421 //adcSumQ_ZNA_lg+= fDigit->GetADCValue(1);
424 GetDigitsData(9)->Fill(fDigit->GetADCValue(0));
425 //GetDigitsData(21)->Fill(fDigit->GetADCValue(1));
428 else if(fDigit->GetSector(0)==5){
429 adcSum_ZPA += fDigit->GetADCValue(0);
430 //adcSum_ZPA_lg += fDigit->GetADCValue(1);
432 if(fDigit->GetSector(1)!=0){
433 adcSumQ_ZPA += fDigit->GetADCValue(0);
434 //adcSumQ_ZPA_lg+= fDigit->GetADCValue(1);
437 GetDigitsData(11)->Fill(fDigit->GetADCValue(0));
438 //GetDigitsData(23)->Fill(fDigit->GetADCValue(1));
443 GetDigitsData(0)->Fill(adcSum_ZNC);
444 GetDigitsData(1)->Fill(adcSum_ZNA);
445 GetDigitsData(2)->Fill(adcSum_ZPC);
446 GetDigitsData(3)->Fill(adcSum_ZPA);
448 GetDigitsData(4)->Fill(adcSumQ_ZNC);
449 GetDigitsData(5)->Fill(adcSumQ_ZNA);
450 GetDigitsData(6)->Fill(adcSumQ_ZPC);
451 GetDigitsData(7)->Fill(adcSumQ_ZPA);
453 /*GetDigitsData(12)->Fill(adcSum_ZNC_lg);
454 GetDigitsData(13)->Fill(adcSum_ZNA_lg);
455 GetDigitsData(14)->Fill(adcSum_ZPC_lg);
456 GetDigitsData(15)->Fill(adcSum_ZPA_lg);
458 GetDigitsData(16)->Fill(adcSumQ_ZNC_lg);
459 GetDigitsData(17)->Fill(adcSumQ_ZNA_lg);
460 GetDigitsData(18)->Fill(adcSumQ_ZPC_lg);
461 GetDigitsData(19)->Fill(adcSumQ_ZPA_lg);*/
464 //____________________________________________________________________________
465 void AliZDCQADataMakerRec::MakeESDs(AliESDEvent * esd)
467 // make QA data from ESDs
470 // Check id histograms already created for this Event Specie
471 if ( ! GetESDsData(0) )
474 AliESDZDC * zdcESD = esd->GetESDZDC();
476 Double32_t * centr_ZNC = zdcESD->GetZNCCentroid();
477 GetESDsData(0)->Fill(centr_ZNC[0], centr_ZNC[1]);
479 Double32_t * centr_ZNA = zdcESD->GetZNACentroid();
480 GetESDsData(1)->Fill(centr_ZNA[0], centr_ZNA[1]);
483 GetESDsData(2)->Fill(esd->GetZDCN1Energy());
484 GetESDsData(3)->Fill(esd->GetZDCN2Energy());
485 GetESDsData(4)->Fill(esd->GetZDCP1Energy());
486 GetESDsData(5)->Fill(esd->GetZDCP2Energy());
488 Double_t sumQZNC=0., sumQZPC=0., sumQZNA=0., sumQZPA=0.;
489 //Double_t sumQZNC_lg=0., sumQZPC_lg=0., sumQZNA_lg=0., sumQZPA_lg=0.;
491 const Double_t *towZNC, *towZPC, *towZNA, *towZPA;
492 //const Double_t *towZNC_lg, *towZPC_lg, *towZNA_lg, *towZPA_lg;
494 towZNC = zdcESD->GetZN1TowerEnergy();
495 towZPC = zdcESD->GetZP1TowerEnergy();
496 towZNA = zdcESD->GetZN2TowerEnergy();
497 towZPA = zdcESD->GetZP2TowerEnergy();
499 /*towZNC_lg = zdcESD->GetZN1TowerEnergyLR();
500 towZPC_lg = zdcESD->GetZP1TowerEnergyLR();
501 towZNA_lg = zdcESD->GetZN2TowerEnergyLR();
502 towZPA_lg = zdcESD->GetZP2TowerEnergyLR();*/
504 for(Int_t i=0; i<5; i++){
506 GetESDsData(10)->Fill(towZNC[i]);
507 GetESDsData(11)->Fill(towZNA[i]);
508 GetESDsData(12)->Fill(towZPC[i]);
509 GetESDsData(13)->Fill(towZPA[i]);
511 /*GetESDsData(22)->Fill(towZNC_lg[i]);
512 GetESDsData(23)->Fill(towZNA_lg[i]);
513 GetESDsData(24)->Fill(towZPC_lg[i]);
514 GetESDsData(25)->Fill(towZPA_lg[i]);*/
517 sumQZNC += towZNC[i];
518 sumQZPC += towZPC[i];
519 sumQZNA += towZNA[i];
520 sumQZPA += towZPA[i];
522 /*sumQZNC_lg += towZNC_lg[i];
523 sumQZPC_lg += towZPC_lg[i];
524 sumQZNA_lg += towZNA_lg[i];
525 sumQZPA_lg += towZPA_lg[i];*/
528 GetESDsData(6)->Fill(sumQZNC);
529 GetESDsData(7)->Fill(sumQZNA);
530 GetESDsData(8)->Fill(sumQZPC);
531 GetESDsData(9)->Fill(sumQZPA);
533 /*GetESDsData(18)->Fill(sumQZNC_lg);
534 GetESDsData(19)->Fill(sumQZNA_lg);
535 GetESDsData(20)->Fill(sumQZPC_lg);
536 GetESDsData(21)->Fill(sumQZPA_lg);*/
539 //____________________________________________________________________________
540 void AliZDCQADataMakerRec::StartOfDetectorCycle()
542 //Detector specific actions at start of cycle
546 //____________________________________________________________________________
547 void AliZDCQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
549 //Detector specific actions at end of cycle
550 // do the QA checking
551 AliQAChecker::Instance()->Run(AliQAv1::kZDC, task, list) ;