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>
31 #include <Riostream.h>
32 // --- Standard library ---
34 // --- AliRoot header files ---
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"
46 ClassImp(AliZDCQADataMakerRec)
48 //____________________________________________________________________________
49 AliZDCQADataMakerRec::AliZDCQADataMakerRec() :
50 AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kZDC), "ZDC Quality Assurance Data Maker"),
56 //____________________________________________________________________________
57 AliZDCQADataMakerRec::AliZDCQADataMakerRec(const AliZDCQADataMakerRec& qadm) :
59 fPedCalibData(qadm.fPedCalibData)
63 SetName((const char*)qadm.GetName());
64 SetTitle((const char*)qadm.GetTitle());
67 //__________________________________________________________________
68 AliZDCQADataMakerRec& AliZDCQADataMakerRec::operator = (const AliZDCQADataMakerRec& qadm )
71 this->~AliZDCQADataMakerRec();
72 new(this) AliZDCQADataMakerRec(qadm);
76 //____________________________________________________________________________
77 AliZDCQADataMakerRec::~AliZDCQADataMakerRec()
79 if(fPedCalibData && !(AliCDBManager::Instance()->GetCacheFlag())){
85 //____________________________________________________________________________
86 AliZDCPedestals* AliZDCQADataMakerRec::GetPedCalibData() const
89 // Retrieving pedestal calibration object from OCDB
90 AliCDBEntry *entry = AliCDBManager::Instance()->Get("ZDC/Calib/Pedestals");
91 if(!entry) AliWarning("No calibration data loaded!");
93 AliZDCPedestals *calibdata = (AliZDCPedestals*) (entry->GetObject());
94 if(!calibdata) AliFatal("Wrong calibration object in calibration file!");
100 //____________________________________________________________________________
101 void AliZDCQADataMakerRec::InitDigits()
103 // create Digits histograms in Digits subdir
105 const Bool_t expert = kTRUE ;
106 const Bool_t image = kTRUE ;
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);
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);
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);
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);
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);
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);
166 //____________________________________________________________________________
167 void AliZDCQADataMakerRec::InitRaws()
169 // create Digits histograms in Digits subdir
170 const Bool_t expert = kTRUE ;
171 const Bool_t image = kTRUE ;
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,7., 1207.);
178 TH1F * hZEM2Spectrum = new TH1F("hZEM2Spectrum","ZEM2 spectrum;Amplitude [ADC counts];Counts",100,7., 1207.);
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);
186 TH2F * hZNCpmCvsPMq = new TH2F("hZNCpmCvsPMq", "ZNC;PMC [ADC counts];Sum(PMQ) [ADC counts]",50,7.,1207.,50,7., 1207.);
187 TH2F * hZPCpmCvsPMq = new TH2F("hZPCpmCvsPMq", "ZPC;PMC [ADC counts];Sum(PMQ) [ADC counts]",50,7.,1207.,50,7., 1207.);
188 TH2F * hZNApmCvsPMq = new TH2F("hZNApmCvsPMq", "ZNA;PMC [ADC counts];Sum(PMQ) [ADC counts]",50,7.,1207.,50,7., 1207.);
189 TH2F * hZPApmCvsPMq = new TH2F("hZPApmCvsPMq", "ZPA;PMC [ADC counts];Sum(PMQ) [ADC counts]",50,7.,1207.,50,7., 1207.);
190 Add2RawsList(hZNCpmCvsPMq, 6, !expert, image);
191 Add2RawsList(hZNApmCvsPMq, 7, !expert, image);
192 Add2RawsList(hZPCpmCvsPMq, 8, !expert, image);
193 Add2RawsList(hZPApmCvsPMq, 9, !expert, image);
195 TH1F * hRawPMCZNC = new TH1F("hRawPMCZNC", "Raw ZNC PMC;Amplitude [ADC counts];Counts",100, 7., 1207.);
196 TH1F * hRawPMCZNA = new TH1F("hRawPMCZNA", "Raw ZNA PMC;Amplitude [ADC counts];Counts",100, 7., 1207.);
197 TH1F * hRawPMCZPC = new TH1F("hRawPMCZPC", "Raw ZPC PMC;Amplitude [ADC counts];Counts",100, 7., 1207.);
198 TH1F * hRawPMCZPA = new TH1F("hRawPMCZPA", "Raw ZPA PMC;Amplitude [ADC counts];Counts",100, 7., 1207.);
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);
212 TH1F * hRawTDCZEM1 = new TH1F("hRawTDCZEM1", "Raw TDC ZEM1;TDC [ns]",200, -100., -50.);
213 Add2RawsList(hRawTDCZEM1, 18, !expert, image);
214 TH1F * hRawTDCZPC = new TH1F("hRawTDCZPC", "Raw TDC ZPC;TDC [ns]",200, -100., -50.);
215 Add2RawsList(hRawTDCZPC, 19, !expert, image);
219 //____________________________________________________________________________
220 void AliZDCQADataMakerRec::InitRecPoints()
222 // create Digits histograms in Digits subdir
223 const Bool_t expert = kTRUE ;
224 const Bool_t image = kTRUE ;
226 // ------------------- HIGH GAIN CHAIN ---------------------------
227 TH1F * hRecZNCTot = new TH1F("hRecZNCTot", "Rec signal in ZNC;Amplitude [ADC counts];Counts", 100, 0., 2000.);
228 TH1F * hRecZNATot = new TH1F("hRecZNATot", "Rec signal in ZNA;Amplitude [ADC counts];Counts", 100, 0., 2000.);
229 TH1F * hRecZPCTot = new TH1F("hRecZPCTot", "Rec signal in ZPC;Amplitude [ADC counts];Counts", 100, 0., 8000.);
230 TH1F * hRecZPATot = new TH1F("hRecZPATot", "Rec signal in ZPA;Amplitude [ADC counts];Counts", 100, 0., 8000.);
231 Add2RecPointsList(hRecZNCTot, 0, !expert, image);
232 Add2RecPointsList(hRecZNATot, 1, !expert, image);
233 Add2RecPointsList(hRecZPCTot, 2, !expert, image);
234 Add2RecPointsList(hRecZPATot, 3, !expert, image);
236 TH1F * hRecSumQZNC = new TH1F("hRecSumQZNC", "Rec summed 4 ZNC quadrants;Amplitude [ADC counts];Counts",100, 0., 2000.);
237 TH1F * hRecSumQZNA = new TH1F("hRecSumQZNA", "Rec summed 4 ZNA quadrants;Amplitude [ADC counts];Counts",100, 0., 2000.);
238 TH1F * hRecSumQZPC = new TH1F("hRecSumQZPC", "Rec summed 4 ZPC quadrants;Amplitude [ADC counts];Counts",100, 0., 2000.);
239 TH1F * hRecSumQZPA = new TH1F("hRecSumQZPA", "Rec summed 4 ZPA quadrants;Amplitude [ADC counts];Counts",100, 0., 2000.);
240 Add2RecPointsList(hRecSumQZNC, 4, expert, !image);
241 Add2RecPointsList(hRecSumQZNA, 5, expert, !image);
242 Add2RecPointsList(hRecSumQZPC, 6, expert, !image);
243 Add2RecPointsList(hRecSumQZPA, 7, expert, !image);
245 TH1F * hRecPMCZNC = new TH1F("hRecPMCZNC", "Rec common ZNC PMT;Amplitude [ADC counts];Counts",100, 0., 2000.);
246 TH1F * hRecPMCZNA = new TH1F("hRecPMCZNA", "Rec common ZNA PMT;Amplitude [ADC counts];Counts",100, 0., 2000.);
247 TH1F * hRecPMCZPC = new TH1F("hRecPMCZPC", "Rec common ZPC PMT;Amplitude [ADC counts];Counts",100, 0., 2000.);
248 TH1F * hRecPMCZPA = new TH1F("hRecPMCZPA", "Rec common ZPA PMT;Amplitude [ADC counts];Counts",100, 0., 2000.);
249 Add2RecPointsList(hRecPMCZNC, 8 , expert, !image);
250 Add2RecPointsList(hRecPMCZNA, 9 , expert, !image);
251 Add2RecPointsList(hRecPMCZPC, 10, expert, !image);
252 Add2RecPointsList(hRecPMCZPA, 11, expert, !image);
255 //____________________________________________________________________________
256 void AliZDCQADataMakerRec::InitESDs()
258 //Booking ESDs histograms
260 const Bool_t expert = kTRUE ;
261 const Bool_t image = kTRUE ;
263 TH2F * hZNCcentr = new TH2F("hZNCcentr", "Centroid in ZNC;X [cm];Y[cm]", 100, -5.,5.,100,-5.,5.);
264 TH2F * hZNAcentr = new TH2F("hZNAcentr", "Centroid in ZNA;X [cm];Y[cm]", 100, -5.,5.,100,-5.,5.);
265 Add2ESDsList(hZNCcentr, 0, !expert, image);
266 Add2ESDsList(hZNAcentr, 1, !expert, image);
268 // ------------------- HIGH GAIN CHAIN ---------------------------
269 TH1F * hESDZNCTot = new TH1F("hESDZNCTot", "Energy in ZNC", 100, 0., 4000.);
270 TH1F * hESDZNATot = new TH1F("hESDZNATot", "Energy in ZNA", 100, 0., 4000.);
271 TH1F * hESDZPCTot = new TH1F("hESDZPCTot", "Energy in ZPC", 100, 0., 4000.);
272 TH1F * hESDZPATot = new TH1F("hESDZPATot", "Energy in ZPA", 100, 0., 4000.);
273 Add2ESDsList(hESDZNCTot, 2, !expert, image);
274 Add2ESDsList(hESDZNATot, 3, !expert, image);
275 Add2ESDsList(hESDZPCTot, 4, !expert, image);
276 Add2ESDsList(hESDZPATot, 5, !expert, image);
278 TH1F * hESDZEM1 = new TH1F("hESDZEM1", "Energy in ZEM1", 100, 0., 2000.);
279 TH1F * hESDZEM2 = new TH1F("hESDZEM2", "Energy in ZEM2", 100, 0., 2000.);
280 Add2ESDsList(hESDZEM1,6, !expert, image);
281 Add2ESDsList(hESDZEM2,7, !expert, image);
283 TH1F * hESDSumQZNC = new TH1F("hESDSumQZNC", "Sum of 4 ZNC energy",100, 0., 2000.);
284 TH1F * hESDSumQZNA = new TH1F("hESDSumQZNA", "Sum of 4 ZNA energy",100, 0., 2000.);
285 TH1F * hESDSumQZPC = new TH1F("hESDSumQZPC", "Sum of 4 ZPC energy",100, 0., 2000.);
286 TH1F * hESDSumQZPA = new TH1F("hESDSumQZPA", "Sum of 4 ZPA energy",100, 0., 2000.);
287 Add2ESDsList(hESDSumQZNC, 8, expert, !image);
288 Add2ESDsList(hESDSumQZNA, 9, expert, !image);
289 Add2ESDsList(hESDSumQZPC, 10, expert, !image);
290 Add2ESDsList(hESDSumQZPA, 11, expert, !image);
292 TH1F * hESDPMCZNC = new TH1F("hESDPMCZNC", "Energy in ZNC PMC",100, 0., 2000.);
293 TH1F * hESDPMCZNA = new TH1F("hESDPMCZNA", "Energy in ZNA PMC",100, 0., 2000.);
294 TH1F * hESDPMCZPC = new TH1F("hESDPMCZPC", "Energy in ZPC PMC",100, 0., 2000.);
295 TH1F * hESDPMCZPA = new TH1F("hESDPMCZPA", "Energy in ZPA PMC",100, 0., 2000.);
296 Add2ESDsList(hESDPMCZNC, 12, expert, !image);
297 Add2ESDsList(hESDPMCZNA, 13, expert, !image);
298 Add2ESDsList(hESDPMCZPC, 14, expert, !image);
299 Add2ESDsList(hESDPMCZPA, 15, expert, !image);
301 // ------------------- LOW GAIN CHAIN ---------------------------
302 TH1F * hESDSumQZNClg = new TH1F("hESDSumQZNClg", "Sum of 4 lg ZNC sectors",100, 0., 4000.);
303 TH1F * hESDSumQZNAlg = new TH1F("hESDSumQZNAlg", "Sum of 4 lg ZNA sectors",100, 0., 4000.);
304 TH1F * hESDSumQZPClg = new TH1F("hESDSumQZPClg", "Sum of 4 lg ZPC sectors",100, 0., 4000.);
305 TH1F * hESDSumQZPAlg = new TH1F("hESDSumQZPAlg", "Sum of 4 lg ZPA sectors",100, 0., 4000.);
306 Add2ESDsList(hESDSumQZNClg, 16, expert, !image);
307 Add2ESDsList(hESDSumQZNAlg, 17, expert, !image);
308 Add2ESDsList(hESDSumQZPClg, 18, expert, !image);
309 Add2ESDsList(hESDSumQZPAlg, 19, expert, !image);
311 TH1F * hESDPMCZNClg = new TH1F("hESDPMCZNClg", "Signal in common ZNC lg PMT",100, 0., 4000.);
312 TH1F * hESDPMCZNAlg = new TH1F("hESDPMCZNAlg", "Signal in common ZNA lg PMT",100, 0., 4000.);
313 TH1F * hESDPMCZPClg = new TH1F("hESDPMCZPClg", "Signal in common ZPC lg PMT",100, 0., 4000.);
314 TH1F * hESDPMCZPAlg = new TH1F("hESDPMCZPAlg", "Signal in common ZPA lg PMT",100, 0., 4000.);
315 Add2ESDsList(hESDPMCZNClg, 20, expert, !image);
316 Add2ESDsList(hESDPMCZNAlg, 21, expert, !image);
317 Add2ESDsList(hESDPMCZPClg, 22, expert, !image);
318 Add2ESDsList(hESDPMCZPAlg, 23, expert, !image);
321 //___________________________________________________________________________
322 void AliZDCQADataMakerRec::MakeDigits(TTree *digitTree)
324 // makes data from Digit Tree
325 if(!GetDigitsData(0)) InitDigits();
328 AliError("Can't get ZDC digit tree!!");
332 TBranch * branch = digitTree->GetBranch("ZDC");
334 AliError("ZDC branch in digit tree not found");
338 AliZDCDigit *digit = 0x0;
339 branch->SetAddress(&digit);
341 Float_t adcSum_ZNC=0., adcSum_ZNA=0., adcSum_ZPC=0., adcSum_ZPA=0.;
342 Float_t adcSumQ_ZNC=0., adcSumQ_ZNA=0., adcSumQ_ZPC=0., adcSumQ_ZPA=0.;
343 Float_t adcSum_ZNC_lg=0., adcSum_ZNA_lg=0., adcSum_ZPC_lg=0., adcSum_ZPA_lg=0.;
344 Float_t adcSumQ_ZNC_lg=0., adcSumQ_ZNA_lg=0., adcSumQ_ZPC_lg=0., adcSumQ_ZPA_lg=0.;
346 Int_t ndig = digitTree->GetEntries();
347 for(Int_t i=0; i<ndig; i++){
350 if(digit->GetSector(0)==1 && digit->GetSector(1)!=5){
351 adcSum_ZNC += digit->GetADCValue(0);
352 adcSum_ZNC_lg += digit->GetADCValue(1);
354 if(digit->GetSector(1)!=0){
355 adcSumQ_ZNC += digit->GetADCValue(0);
356 adcSumQ_ZNC_lg+= digit->GetADCValue(1);
359 GetDigitsData(8)->Fill(digit->GetADCValue(0));
360 GetDigitsData(20)->Fill(digit->GetADCValue(1));
363 else if(digit->GetSector(0)==2){
364 adcSum_ZPC += digit->GetADCValue(0);
365 adcSum_ZPC_lg += digit->GetADCValue(1);
367 if(digit->GetSector(1)!=0){
368 adcSumQ_ZPC += digit->GetADCValue(0);
369 adcSumQ_ZPC_lg+= digit->GetADCValue(1);
372 GetDigitsData(10)->Fill(digit->GetADCValue(0));
373 GetDigitsData(22)->Fill(digit->GetADCValue(1));
376 else if(digit->GetSector(0)==4 && digit->GetSector(1)!=5){
377 adcSum_ZNA += digit->GetADCValue(0);
378 adcSum_ZNA_lg += digit->GetADCValue(1);
380 if(digit->GetSector(1)!=0){
381 adcSumQ_ZNA += digit->GetADCValue(0);
382 adcSumQ_ZNA_lg+= digit->GetADCValue(1);
385 GetDigitsData(9)->Fill(digit->GetADCValue(0));
386 GetDigitsData(21)->Fill(digit->GetADCValue(1));
389 else if(digit->GetSector(0)==5){
390 adcSum_ZPA += digit->GetADCValue(0);
391 adcSum_ZPA_lg += digit->GetADCValue(1);
393 if(digit->GetSector(1)!=0){
394 adcSumQ_ZPA += digit->GetADCValue(0);
395 adcSumQ_ZPA_lg+= digit->GetADCValue(1);
398 GetDigitsData(11)->Fill(digit->GetADCValue(0));
399 GetDigitsData(23)->Fill(digit->GetADCValue(1));
404 GetDigitsData(0)->Fill(adcSum_ZNC);
405 GetDigitsData(1)->Fill(adcSum_ZNA);
406 GetDigitsData(2)->Fill(adcSum_ZPC);
407 GetDigitsData(3)->Fill(adcSum_ZPA);
409 GetDigitsData(4)->Fill(adcSumQ_ZNC);
410 GetDigitsData(5)->Fill(adcSumQ_ZNA);
411 GetDigitsData(6)->Fill(adcSumQ_ZPC);
412 GetDigitsData(7)->Fill(adcSumQ_ZPA);
414 GetDigitsData(12)->Fill(adcSum_ZNC_lg);
415 GetDigitsData(13)->Fill(adcSum_ZNA_lg);
416 GetDigitsData(14)->Fill(adcSum_ZPC_lg);
417 GetDigitsData(15)->Fill(adcSum_ZPA_lg);
419 GetDigitsData(16)->Fill(adcSumQ_ZNC_lg);
420 GetDigitsData(17)->Fill(adcSumQ_ZNA_lg);
421 GetDigitsData(18)->Fill(adcSumQ_ZPC_lg);
422 GetDigitsData(19)->Fill(adcSumQ_ZPA_lg);
430 //____________________________________________________________________________
431 void AliZDCQADataMakerRec::MakeRaws(AliRawReader *rawReader)
433 // Filling Raws QA histos
435 // Checking the event type
436 // if (rawReader->GetType()!=7){
438 // Check if histograms already created for this Event Specie
439 if(!GetRawsData(0)) InitRaws();
441 // Parameters for mean value pedestal subtraction
443 Float_t meanPed[2*kNch];
444 for(Int_t jj=0; jj<2*kNch; jj++) meanPed[jj] = fPedCalibData->GetMeanPed(jj);
446 Float_t zncSignal=0., znaSignal=0., zpcSignal=0., zpaSignal=0.;
447 Float_t zncSumQ=0., znaSumQ=0., zpcSumQ=0., zpaSumQ=0.;
448 Float_t zncpmC=0., znapmC=0., zpcpmC=0., zpapmC=0.;
449 Bool_t isZNCFired=kFALSE, isZPCFired=kFALSE, isZNAFired=kFALSE, isZPAFired=kFALSE;
450 Int_t indZNC=0, indZNA=0, indZPC=0, indZPA=0;
451 Float_t zpcTDC[10], zem1TDC[10];
452 for(Int_t i=0; i<10; i++){
453 zpcTDC[i] = zem1TDC[i] = -999.;
456 Int_t iMultZPCTDC=0, iMultZEMTDC=0;
459 AliZDCRawStream stream(rawReader);
460 while(stream.Next()){
462 if(stream.IsADCDataWord() &&
463 (stream.GetADCModule()==0 || stream.GetADCModule()==1)){
465 Int_t det = stream.GetSector(0);
466 Int_t quad = stream.GetSector(1);
467 Int_t gain = stream.GetADCGain();
470 // Stuff for pedestal subtraction
471 if(quad != 5){ // ZDCs (not reference PTMs)
475 Float_t pedSubVal = (Float_t) (stream.GetADCValue()-meanPed[pedindex]);
476 zncSignal += pedSubVal;
478 if(quad!=0) zncSumQ += pedSubVal;
481 GetRawsData(10)->Fill(zncpmC);
489 Float_t pedSubVal = (Float_t) (stream.GetADCValue()-meanPed[pedindex]);
490 zpcSignal += pedSubVal;
492 if(quad!=0) zpcSumQ += pedSubVal;
495 GetRawsData(12)->Fill(zpcpmC);
504 Float_t pedSubVal = (Float_t) (stream.GetADCValue()-meanPed[pedindex]);
505 GetRawsData(4)->Fill(pedSubVal);
510 Float_t pedSubVal = (Float_t) (stream.GetADCValue()-meanPed[pedindex]);
511 GetRawsData(5)->Fill(pedSubVal);
518 Float_t pedSubVal = (Float_t) (stream.GetADCValue()-meanPed[pedindex]);
519 znaSignal += pedSubVal;
521 if(quad!=0) znaSumQ += pedSubVal;
524 GetRawsData(11)->Fill(znapmC);
532 Float_t pedSubVal = (Float_t) (stream.GetADCValue()-meanPed[pedindex]);
533 zpaSignal += pedSubVal;
535 if(quad!=0) zpaSumQ += pedSubVal;
538 GetRawsData(13)->Fill(zpapmC);
545 if(isZNCFired && indZNC==5){
546 GetRawsData(0)->Fill(zncSignal);
547 GetRawsData(6)->Fill(zncSumQ, zncpmC);
548 GetRawsData(14)->Fill(zncSumQ);
550 if(isZPCFired && indZPC==5){
551 GetRawsData(2)->Fill(zpcSignal);
552 GetRawsData(8)->Fill(zpcSumQ, zpcpmC);
553 GetRawsData(16)->Fill(zpcSumQ);
555 if(isZNAFired && indZNA==5){
556 GetRawsData(1)->Fill(znaSignal);
557 GetRawsData(7)->Fill(znaSumQ, znapmC);
558 GetRawsData(15)->Fill(znaSumQ);
560 if(isZPAFired && indZPA==5){
561 GetRawsData(3)->Fill(zpaSignal);
562 GetRawsData(9)->Fill(zpaSumQ, zpapmC);
563 GetRawsData(17)->Fill(zpaSumQ);
567 zncSignal = zncSumQ = zncpmC = 0;
568 isZNCFired=kFALSE; indZNC=0;
571 zpcSignal = zpcSumQ = zpcpmC = 0;
572 isZPCFired=kFALSE; indZPC=0;
575 znaSignal = znaSumQ = znapmC = 0;
576 isZNAFired=kFALSE; indZNA=0;
579 zpaSignal = zpaSumQ = zpapmC = 0;
580 isZPAFired=kFALSE; indZPA=0;
583 } //IsADCDataWord && signal ADCs
584 else if(stream.IsZDCTDCDatum()){
585 if(stream.GetChannel()==3){
586 zpcTDC[iMultZPCTDC] = (0.025*(stream.GetZDCTDCDatum()));
589 else if(stream.GetChannel()==8){
590 zem1TDC[iMultZEMTDC] = (0.025*(stream.GetZDCTDCDatum()));
593 else if(stream.GetChannel()==16) tdcL0 = (0.025*(stream.GetZDCTDCDatum()));
595 if(stream.GetChannel()==16 && tdcL0!=-999.){
596 for(Int_t iHit=0; iHit<10; iHit++){
597 if(zpcTDC[iHit]!=-999.){
598 Float_t diffZPC = zpcTDC[iHit]-tdcL0;
599 GetRawsData(19)->Fill(diffZPC);
601 if(zem1TDC[iHit]!=-999.){
602 Float_t diffZEM1 = zem1TDC[iHit]-tdcL0;
603 GetRawsData(18)->Fill(diffZEM1);
608 for(Int_t i=0; i<10; i++){
609 zpcTDC[i] = zem1TDC[i] = -999.;
615 // } // check on event type
617 // AliDebug(1,Form("Skipping non-physics event for QA -> event type %d \n", rawReader->GetType()));
621 //____________________________________________________________________________
622 void AliZDCQADataMakerRec::MakeRecPoints(TTree * clustersTree)
624 // Filling QA histos from RecPoints
626 TBranch *branch = clustersTree->GetBranch("ZDC");
628 AliError("Can't get the ZDC branch for rec points!");
632 if(!GetRecPointsData(0)) InitRecPoints() ;
634 Float_t sum_ZNC=0., sum_ZNA=0., sum_ZPC=0., sum_ZPA=0.;
635 Float_t sumQ_ZNC=0., sumQ_ZNA=0., sumQ_ZPC=0., sumQ_ZPA=0.;
638 AliZDCReco* preco = &reco;
639 clustersTree->SetBranchAddress("ZDC", &preco);
641 clustersTree->GetEntry(0);
642 for(Int_t i=0; i<5; i++){
643 sum_ZNC += reco.GetZN1HREnTow(i);
644 sum_ZPC += reco.GetZN2HREnTow(i);
645 sum_ZNA += reco.GetZP1HREnTow(i);
646 sum_ZPA += reco.GetZP2HREnTow(i);
648 GetRecPointsData(8)->Fill(reco.GetZN1HREnTow(i));
649 GetRecPointsData(9)->Fill(reco.GetZN2HREnTow(i));
650 GetRecPointsData(10)->Fill(reco.GetZP1HREnTow(i));
651 GetRecPointsData(11)->Fill(reco.GetZP2HREnTow(i));
654 sumQ_ZNC += reco.GetZN1HREnTow(i);
655 sumQ_ZPC += reco.GetZN2HREnTow(i);
656 sumQ_ZNA += reco.GetZP1HREnTow(i);
657 sumQ_ZPA += reco.GetZP2HREnTow(i);
661 GetRecPointsData(0)->Fill(sum_ZNC);
662 GetRecPointsData(1)->Fill(sum_ZNA);
663 GetRecPointsData(2)->Fill(sum_ZPC);
664 GetRecPointsData(3)->Fill(sum_ZPA);
666 GetRecPointsData(4)->Fill(sumQ_ZNC);
667 GetRecPointsData(5)->Fill(sumQ_ZNA);
668 GetRecPointsData(6)->Fill(sumQ_ZPC);
669 GetRecPointsData(7)->Fill(sumQ_ZPA);
673 //____________________________________________________________________________
674 void AliZDCQADataMakerRec::MakeESDs(AliESDEvent * esd)
676 // make QA data from ESDs
679 // Check id histograms already created for this Event Specie
680 if(!GetESDsData(0)) InitESDs() ;
682 AliESDZDC * zdcESD = esd->GetESDZDC();
684 TString beamType = esd->GetBeamType();
685 Double_t centr_ZNC[2]={999.,999}, centr_ZNA[2]={999.,999};
686 if(((beamType.CompareTo("pp"))==0) || ((beamType.CompareTo("p-p"))==0)
687 ||((beamType.CompareTo("PP"))==0) || ((beamType.CompareTo("P-P"))==0)){
688 zdcESD->GetZNCentroidInpp(centr_ZNC, centr_ZNA);
690 else if((beamType.CompareTo("A-A")) == 0 || (beamType.CompareTo("Pb-Pb")) == 0){
691 Float_t beamEne = esd->GetBeamEnergy();
692 zdcESD->GetZNCentroidInPbPb(beamEne, centr_ZNC, centr_ZNA);
694 else printf("\n WARNING!!! AliZDCQADataMakerRec::MakeESDs: can't calculate centroids for beam type: %s\n\n",beamType.Data());
695 GetESDsData(0)->Fill(centr_ZNC[0], centr_ZNC[1]);
696 GetESDsData(1)->Fill(centr_ZNA[0], centr_ZNA[1]);
699 GetESDsData(2)->Fill(esd->GetZDCN1Energy());
700 GetESDsData(3)->Fill(esd->GetZDCN2Energy());
701 GetESDsData(4)->Fill(esd->GetZDCP1Energy());
702 GetESDsData(5)->Fill(esd->GetZDCP2Energy());
703 GetESDsData(6)->Fill(esd->GetZDCEMEnergy(0));
704 GetESDsData(7)->Fill(esd->GetZDCEMEnergy(1));
706 Double_t sumQZNC=0., sumQZPC=0., sumQZNA=0., sumQZPA=0.;
707 Double_t sumQZNC_lg=0., sumQZPC_lg=0., sumQZNA_lg=0., sumQZPA_lg=0.;
709 const Double_t *towZNC, *towZPC, *towZNA, *towZPA;
710 const Double_t *towZNC_lg, *towZPC_lg, *towZNA_lg, *towZPA_lg;
712 towZNC = zdcESD->GetZN1TowerEnergy();
713 towZPC = zdcESD->GetZP1TowerEnergy();
714 towZNA = zdcESD->GetZN2TowerEnergy();
715 towZPA = zdcESD->GetZP2TowerEnergy();
717 towZNC_lg = zdcESD->GetZN1TowerEnergyLR();
718 towZPC_lg = zdcESD->GetZP1TowerEnergyLR();
719 towZNA_lg = zdcESD->GetZN2TowerEnergyLR();
720 towZPA_lg = zdcESD->GetZP2TowerEnergyLR();
722 for(Int_t i=0; i<5; i++){
724 GetESDsData(12)->Fill(towZNC[i]);
725 GetESDsData(13)->Fill(towZNA[i]);
726 GetESDsData(14)->Fill(towZPC[i]);
727 GetESDsData(15)->Fill(towZPA[i]);
729 GetESDsData(20)->Fill(towZNC_lg[i]);
730 GetESDsData(21)->Fill(towZNA_lg[i]);
731 GetESDsData(22)->Fill(towZPC_lg[i]);
732 GetESDsData(23)->Fill(towZPA_lg[i]);
735 sumQZNC += towZNC[i];
736 sumQZPC += towZPC[i];
737 sumQZNA += towZNA[i];
738 sumQZPA += towZPA[i];
740 sumQZNC_lg += towZNC_lg[i];
741 sumQZPC_lg += towZPC_lg[i];
742 sumQZNA_lg += towZNA_lg[i];
743 sumQZPA_lg += towZPA_lg[i];
746 GetESDsData(8)->Fill(sumQZNC);
747 GetESDsData(9)->Fill(sumQZNA);
748 GetESDsData(10)->Fill(sumQZPC);
749 GetESDsData(11)->Fill(sumQZPA);
751 GetESDsData(16)->Fill(sumQZNC_lg);
752 GetESDsData(17)->Fill(sumQZNA_lg);
753 GetESDsData(18)->Fill(sumQZPC_lg);
754 GetESDsData(19)->Fill(sumQZPA_lg);
757 //____________________________________________________________________________
758 void AliZDCQADataMakerRec::StartOfDetectorCycle()
760 //Detector specific actions at start of cycle
762 fPedCalibData = GetPedCalibData();
766 //____________________________________________________________________________
767 void AliZDCQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
769 //Detector specific actions at end of cycle
770 // do the QA checking
771 if( task == AliQAv1::kRAWS){
772 if (!GetRawsData(4) || !GetRawsData(5) || !GetRawsData(6) || !GetRawsData(7) ||
773 !GetRawsData(8) || !GetRawsData(9) || !GetRawsData(10) || !GetRawsData(11) ||
774 !GetRawsData(12) || !GetRawsData(13)) {
775 printf(" WARNING!!! AliZDCQADataMaker Rec -> No histogram for DQM found!\n");
779 GetRawsData(6)->SetOption("colz");
780 GetRawsData(7)->SetOption("colz");
781 GetRawsData(8)->SetOption("colz");
782 GetRawsData(9)->SetOption("colz");
784 GetRawsData(4)->SetLineColor(kAzure+1); GetRawsData(4)->SetLineWidth(2);
785 GetRawsData(5)->SetLineColor(kAzure+2); GetRawsData(5)->SetLineWidth(2);
786 GetRawsData(10)->SetLineColor(kBlue); GetRawsData(10)->SetLineWidth(2);
787 GetRawsData(11)->SetLineColor(kBlue+1); GetRawsData(11)->SetLineWidth(2);
788 GetRawsData(12)->SetLineColor(kBlue+2); GetRawsData(12)->SetLineWidth(2);
789 GetRawsData(13)->SetLineColor(kBlue+3); GetRawsData(13)->SetLineWidth(2);
791 /*if(((GetRawsData(4))->GetEntries())>0) GetRawsData(4)->SetDrawOption("LOGY");
792 if(((GetRawsData(5))->GetEntries())>0) GetRawsData(5)->SetDrawOption("LOGY");
793 if(((GetRawsData(10))->GetEntries())>0) GetRawsData(10)->SetDrawOption("LOGY");
794 if(((GetRawsData(11))->GetEntries())>0) GetRawsData(11)->SetDrawOption("LOGY");
795 if(((GetRawsData(12))->GetEntries())>0) GetRawsData(12)->SetDrawOption("LOGY");
796 if(((GetRawsData(13))->GetEntries())>0) GetRawsData(13)->SetDrawOption("LOGY");*/
800 AliQAChecker::Instance()->Run(AliQAv1::kZDC, task, list) ;