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 "AliZDCQADataMakerRec.h"
37 #include "AliZDCRawStream.h"
38 #include "AliESDZDC.h"
39 #include "AliESDEvent.h"
41 ClassImp(AliZDCQADataMakerRec)
43 //____________________________________________________________________________
44 AliZDCQADataMakerRec::AliZDCQADataMakerRec() :
45 AliQADataMakerRec(AliQA::GetDetName(AliQA::kZDC), "ZDC Quality Assurance Data Maker")
50 //____________________________________________________________________________
51 AliZDCQADataMakerRec::AliZDCQADataMakerRec(const AliZDCQADataMakerRec& qadm) :
55 SetName((const char*)qadm.GetName());
56 SetTitle((const char*)qadm.GetTitle());
59 //__________________________________________________________________
60 AliZDCQADataMakerRec& AliZDCQADataMakerRec::operator = (const AliZDCQADataMakerRec& qadm )
63 this->~AliZDCQADataMakerRec();
64 new(this) AliZDCQADataMakerRec(qadm);
68 //____________________________________________________________________________
70 void AliZDCQADataMakerRec::InitRaws()
72 // create Digits histograms in Digits subdir
74 // ------------------- HIGH GAIN CHAIN ---------------------------
75 TH1F * hRawZNCTot = new TH1F("hRawZNCTot", "Raw signal in ZNC", 100, 0., 6000.);
76 TH1F * hRawZNATot = new TH1F("hRawZNATot", "Raw signal in ZNA", 100, 0., 6000.);
77 TH1F * hRawZPCTot = new TH1F("hRawZPCTot", "Raw signal in ZPC", 100, 0., 10000.);
78 TH1F * hRawZPATot = new TH1F("hRawZPATot", "Raw signal in ZPA", 100, 0., 10000.);
79 Add2RawsList(hRawZNCTot, 0);
80 Add2RawsList(hRawZNATot, 1);
81 Add2RawsList(hRawZPCTot, 2);
82 Add2RawsList(hRawZPATot, 3);
84 TH1F * hRawSumQZNC = new TH1F("hRawSumQZNC", "Raw summed 4 ZNC quadrants",100, 0., 4000.);
85 TH1F * hRawSumQZNA = new TH1F("hRawSumQZNA", "Raw summed 4 ZNA quadrants",100, 0., 4000.);
86 TH1F * hRawSumQZPC = new TH1F("hRawSumQZPC", "Raw summed 4 ZPC quadrants",100, 0., 4000.);
87 TH1F * hRawSumQZPA = new TH1F("hRawSumQZPA", "Raw summed 4 ZPA quadrants",100, 0., 4000.);
88 Add2RawsList(hRawSumQZNC, 4);
89 Add2RawsList(hRawSumQZNA, 5);
90 Add2RawsList(hRawSumQZPC, 6);
91 Add2RawsList(hRawSumQZPA, 7);
93 TH1F * hRawPMCZNC = new TH1F("hRawPMCZNC", "Raw common ZNC PMT",100, 0., 4000.);
94 TH1F * hRawPMCZNA = new TH1F("hRawPMCZNA", "Raw common ZNA PMT",100, 0., 4000.);
95 TH1F * hRawPMCZPC = new TH1F("hRawPMCZPC", "Raw common ZPC PMT",100, 0., 4000.);
96 TH1F * hRawPMCZPA = new TH1F("hRawPMCZPA", "Raw common ZPA PMT",100, 0., 4000.);
97 Add2RawsList(hRawPMCZNC, 8);
98 Add2RawsList(hRawPMCZNA, 9);
99 Add2RawsList(hRawPMCZPC, 10);
100 Add2RawsList(hRawPMCZPA, 11);
102 /* // ------------------- LOW GAIN CHAIN ---------------------------
103 TH1F * hRawZNCTotlg = new TH1F("hRawZNCTotlg", "Rawit lg signal in ZNC", 100, 0., 6000.);
104 TH1F * hRawZNATotlg = new TH1F("hRawZNATotlg", "Rawit lg signal in ZNA", 100, 0., 6000.);
105 TH1F * hRawZPCTotlg = new TH1F("hRawZPCTotlg", "Rawit lg signal in ZPC", 100, 0., 10000.);
106 TH1F * hRawZPATotlg = new TH1F("hRawZPATotlg", "Rawit lg signal in ZPA", 100, 0., 10000.);
107 Add2RawsList(hRawZNCTotlg, 12);
108 Add2RawsList(hRawZNATotlg, 13);
109 Add2RawsList(hRawZPCTotlg, 14);
110 Add2RawsList(hRawZPATotlg, 15);
112 TH1F * hRawSumQZNClg = new TH1F("hRawSumQZNClg", "Raw summed 4 lg ZNC quadrants",100, 0., 4000.);
113 TH1F * hRawSumQZNAlg = new TH1F("hRawSumQZNAlg", "Raw summed 4 lg ZNA quadrants",100, 0., 4000.);
114 TH1F * hRawSumQZPClg = new TH1F("hRawSumQZPClg", "Raw summed 4 lg ZPC quadrants",100, 0., 4000.);
115 TH1F * hRawSumQZPAlg = new TH1F("hRawSumQZPAlg", "Raw summed 4 lg ZPA quadrants",100, 0., 4000.);
116 Add2RawsList(hRawSumQZNClg, 16);
117 Add2RawsList(hRawSumQZNAlg, 17);
118 Add2RawsList(hRawSumQZPClg, 18);
119 Add2RawsList(hRawSumQZPAlg, 19);
121 TH1F * hRawPMCZNClg = new TH1F("hRawPMCZNClg", "Raw common lg ZNC PMT",100, 0., 4000.);
122 TH1F * hRawPMCZNAlg = new TH1F("hRawPMCZNAlg", "Raw common lg ZNA PMT",100, 0., 4000.);
123 TH1F * hRawPMCZPClg = new TH1F("hRawPMCZPClg", "Raw common lg ZPC PMT",100, 0., 4000.);
124 TH1F * hRawPMCZPAlg = new TH1F("hRawPMCZPAlg", "Raw common lg ZPA PMT",100, 0., 4000.);
125 Add2RawsList(hRawPMCZNClg, 20);
126 Add2RawsList(hRawPMCZNAlg, 21);
127 Add2RawsList(hRawPMCZPClg, 22);
128 Add2RawsList(hRawPMCZPAlg, 23);*/
131 //____________________________________________________________________________
132 void AliZDCQADataMakerRec::InitESDs()
134 //Booking ESDs histograms
136 TH2F * hZNC = new TH2F("hZNC", "Centroid in ZNC", 100, -5.,5.,100,-5.,5.);
137 TH2F * hZNA = new TH2F("hZNA", "Centroid in ZNA", 100, -5.,5.,100,-5.,5.);
138 Add2ESDsList(hZNC, 0);
139 Add2ESDsList(hZNA, 1);
141 // ------------------- HIGH GAIN CHAIN ---------------------------
142 TH1F * hESDZNCTot = new TH1F("hESDZNCTot", "Energy in ZNC", 100, 0., 6000.);
143 TH1F * hESDZNATot = new TH1F("hESDZNATot", "Energy in ZNA", 100, 0., 6000.);
144 TH1F * hESDZPCTot = new TH1F("hESDZPCTot", "Energy in ZPC", 100, 0., 10000.);
145 TH1F * hESDZPATot = new TH1F("hESDZPATot", "Energy in ZPA", 100, 0., 10000.);
146 Add2ESDsList(hESDZNCTot, 2);
147 Add2ESDsList(hESDZNATot, 3);
148 Add2ESDsList(hESDZPCTot, 4);
149 Add2ESDsList(hESDZPATot, 5);
151 TH1F * hESDSumQZNC = new TH1F("hESDSumQZNC", "Sum of 4 ZNC energy",100, 0., 4000.);
152 TH1F * hESDSumQZNA = new TH1F("hESDSumQZNA", "Sum of 4 ZNA energy",100, 0., 4000.);
153 TH1F * hESDSumQZPC = new TH1F("hESDSumQZPC", "Sum of 4 ZPC energy",100, 0., 4000.);
154 TH1F * hESDSumQZPA = new TH1F("hESDSumQZPA", "Sum of 4 ZPA energy",100, 0., 4000.);
155 Add2ESDsList(hESDSumQZNC, 6);
156 Add2ESDsList(hESDSumQZNA, 7);
157 Add2ESDsList(hESDSumQZPC, 8);
158 Add2ESDsList(hESDSumQZPA, 9);
160 TH1F * hESDPMCZNC = new TH1F("hESDPMCZNC", "Energy in common ZNC PMT",100, 0., 4000.);
161 TH1F * hESDPMCZNA = new TH1F("hESDPMCZNA", "Energy in common ZNA PMT",100, 0., 4000.);
162 TH1F * hESDPMCZPC = new TH1F("hESDPMCZPC", "Energy in common ZPC PMT",100, 0., 4000.);
163 TH1F * hESDPMCZPA = new TH1F("hESDPMCZPA", "Energy in common ZPA PMT",100, 0., 4000.);
164 Add2ESDsList(hESDPMCZNC, 10);
165 Add2ESDsList(hESDPMCZNA, 11);
166 Add2ESDsList(hESDPMCZPC, 12);
167 Add2ESDsList(hESDPMCZPA, 13);
169 /* // ------------------- LOW GAIN CHAIN ---------------------------
170 TH1F * hESDZNCTotlg = new TH1F("hESDZNCTotlg", "ESD lg signal in ZNC", 100, 0., 6000.);
171 TH1F * hESDZNATotlg = new TH1F("hESDZNATotlg", "ESD lg signal in ZNA", 100, 0., 6000.);
172 TH1F * hESDZPCTotlg = new TH1F("hESDZPCTotlg", "ESD lg signal in ZPC", 100, 0., 10000.);
173 TH1F * hESDZPATotlg = new TH1F("hESDZPATotlg", "ESD lg signal in ZPA", 100, 0., 10000.);
174 Add2ESDsList(hESDZNCTotlg, 14);
175 Add2ESDsList(hESDZNATotlg, 15);
176 Add2ESDsList(hESDZPCTotlg, 16);
177 Add2ESDsList(hESDZPATotlg, 17);
179 TH1F * hESDSumQZNClg = new TH1F("hESDSumQZNClg", "Sum of 4 lg ZNC sectors",100, 0., 4000.);
180 TH1F * hESDSumQZNAlg = new TH1F("hESDSumQZNAlg", "Sum of 4 lg ZNA sectors",100, 0., 4000.);
181 TH1F * hESDSumQZPClg = new TH1F("hESDSumQZPClg", "Sum of 4 lg ZPC sectors",100, 0., 4000.);
182 TH1F * hESDSumQZPAlg = new TH1F("hESDSumQZPAlg", "Sum of 4 lg ZPA sectors",100, 0., 4000.);
183 Add2ESDsList(hESDSumQZNClg, 18);
184 Add2ESDsList(hESDSumQZNAlg, 19);
185 Add2ESDsList(hESDSumQZPClg, 20);
186 Add2ESDsList(hESDSumQZPAlg, 21);
188 TH1F * hESDPMCZNClg = new TH1F("hESDPMCZNClg", "Signal in common ZNC lg PMT",100, 0., 4000.);
189 TH1F * hESDPMCZNAlg = new TH1F("hESDPMCZNAlg", "Signal in common ZNA lg PMT",100, 0., 4000.);
190 TH1F * hESDPMCZPClg = new TH1F("hESDPMCZPClg", "Signal in common ZPC lg PMT",100, 0., 4000.);
191 TH1F * hESDPMCZPAlg = new TH1F("hESDPMCZPAlg", "Signal in common ZPA lg PMT",100, 0., 4000.);
192 Add2ESDsList(hESDPMCZNClg, 22);
193 Add2ESDsList(hESDPMCZNAlg, 23);
194 Add2ESDsList(hESDPMCZPClg, 24);
195 Add2ESDsList(hESDPMCZPAlg, 25);*/
198 //____________________________________________________________________________
200 void AliZDCQADataMakerRec::MakeRaws(AliRawReader *rawReader)
202 // Filling Raws QA histos
204 Float_t sum_ZNC=0., sum_ZNA=0., sum_ZPC=0., sum_ZPA=0.;
205 Float_t sumQ_ZNC=0., sumQ_ZNA=0., sumQ_ZPC=0., sumQ_ZPA=0.;
206 //Float_t sum_ZNC_lg=0., sum_ZNA_lg=0., sum_ZPC_lg=0., sum_ZPA_lg=0.;
207 //Float_t sumQ_ZNC_lg=0., sumQ_ZNA_lg=0., sumQ_ZPC_lg=0., sumQ_ZPA_lg=0.;
209 AliZDCRawStream stream(rawReader);
210 while(stream.Next()){
211 if(stream.IsADCDataWord() &&
212 (stream.GetADCModule()==0 || stream.GetADCModule()==1)){
213 if(stream.GetSector(0)==1){
214 if(stream.GetADCGain()==0){
215 sum_ZNC += stream.GetADCValue();
216 if(stream.GetSector(1)!=0) sumQ_ZNC += stream.GetADCValue();
217 else GetRawsData(8)->Fill(stream.GetADCValue());
220 sum_ZNC_lg += stream.GetADCValue();
221 if(stream.GetSector(1)!=0) sumQ_ZNC_lg += stream.GetADCValue();
222 else GetRawsData(20)->Fill(stream.GetADCValue());
225 else if(stream.GetSector(0)==2){
226 if(stream.GetADCGain()==0){
227 sum_ZPC += stream.GetADCValue();
228 if(stream.GetSector(1)!=0) sumQ_ZPC += stream.GetADCValue();
229 else GetRawsData(10)->Fill(stream.GetADCValue());
232 sum_ZPC_lg += stream.GetADCValue();
233 if(stream.GetSector(1)!=0) sumQ_ZPC_lg += stream.GetADCValue();
234 else GetRawsData(22)->Fill(stream.GetADCValue());
237 else if(stream.GetSector(0)==4){
238 if(stream.GetADCGain()==0){
239 sum_ZNA += stream.GetADCValue();
240 if(stream.GetSector(1)!=0) sumQ_ZNA += stream.GetADCValue();
241 else GetRawsData(9)->Fill(stream.GetADCValue());
244 sum_ZNA_lg += stream.GetADCValue();
245 if(stream.GetSector(1)!=0) sumQ_ZNA_lg += stream.GetADCValue();
246 else GetRawsData(21)->Fill(stream.GetADCValue());
249 else if(stream.GetSector(0)==5){
250 if(stream.GetADCGain()==0){
251 sum_ZPA += stream.GetADCValue();
252 if(stream.GetSector(1)!=0) sumQ_ZPA += stream.GetADCValue();
253 else GetRawsData(11)->Fill(stream.GetADCValue());
256 sum_ZPA_lg += stream.GetADCValue();
257 if(stream.GetSector(1)!=0) sumQ_ZPA_lg += stream.GetADCValue();
258 else GetRawsData(23)->Fill(stream.GetADCValue());
264 GetRawsData(0)->Fill(sum_ZNC);
265 GetRawsData(1)->Fill(sum_ZNA);
266 GetRawsData(2)->Fill(sum_ZPC);
267 GetRawsData(3)->Fill(sum_ZPA);
269 GetRawsData(4)->Fill(sumQ_ZNC);
270 GetRawsData(5)->Fill(sumQ_ZNA);
271 GetRawsData(6)->Fill(sumQ_ZPC);
272 GetRawsData(7)->Fill(sumQ_ZPA);
274 /*GetRawsData(12)->Fill(sum_ZNC_lg);
275 GetRawsData(13)->Fill(sum_ZNA_lg);
276 GetRawsData(14)->Fill(sum_ZPC_lg);
277 GetRawsData(15)->Fill(sum_ZPA_lg);
279 GetRawsData(16)->Fill(sumQ_ZNC_lg);
280 GetRawsData(17)->Fill(sumQ_ZNA_lg);
281 GetRawsData(18)->Fill(sumQ_ZPC_lg);
282 GetRawsData(19)->Fill(sumQ_ZPA_lg);*/
287 //____________________________________________________________________________
288 void AliZDCQADataMakerRec::MakeESDs(AliESDEvent * esd)
290 // make QA data from ESDs
292 AliESDZDC * zdcESD = esd->GetESDZDC();
294 Double32_t * centr_ZNC, * centr_ZNA;
295 Int_t NSpecnC = (Int_t) (esd->GetZDCN1Energy()/2.7);
297 centr_ZNC = zdcESD->GetZNCCentroid(NSpecnC);
298 GetESDsData(0)->Fill(centr_ZNC[0], centr_ZNC[1]);
300 Int_t NSpecnA = (Int_t) (esd->GetZDCN2Energy()/2.7);
302 centr_ZNA = zdcESD->GetZNACentroid(NSpecnA);
303 GetESDsData(1)->Fill(centr_ZNA[0], centr_ZNA[1]);
306 GetESDsData(2)->Fill(esd->GetZDCN1Energy());
307 GetESDsData(3)->Fill(esd->GetZDCN2Energy());
308 GetESDsData(4)->Fill(esd->GetZDCP1Energy());
309 GetESDsData(5)->Fill(esd->GetZDCP2Energy());
311 Double_t sumQZNC=0., sumQZPC=0., sumQZNA=0., sumQZPA=0.;
312 //Double_t sumQZNC_lg=0., sumQZPC_lg=0., sumQZNA_lg=0., sumQZPA_lg=0.;
314 const Double_t *towZNC, *towZPC, *towZNA, *towZPA;
315 //const Double_t *towZNC_lg, *towZPC_lg, *towZNA_lg, *towZPA_lg;
317 towZNC = zdcESD->GetZN1TowerEnergy();
318 towZPC = zdcESD->GetZP1TowerEnergy();
319 towZNA = zdcESD->GetZN2TowerEnergy();
320 towZPA = zdcESD->GetZP2TowerEnergy();
322 /*towZNC_lg = zdcESD->GetZN1TowerEnergyLR();
323 towZPC_lg = zdcESD->GetZP1TowerEnergyLR();
324 towZNA_lg = zdcESD->GetZN2TowerEnergyLR();
325 towZPA_lg = zdcESD->GetZP2TowerEnergyLR();*/
327 for(Int_t i=0; i<5; i++){
329 GetESDsData(10)->Fill(towZNC[i]);
330 GetESDsData(11)->Fill(towZNA[i]);
331 GetESDsData(12)->Fill(towZPC[i]);
332 GetESDsData(13)->Fill(towZPA[i]);
334 /*GetESDsData(22)->Fill(towZNC_lg[i]);
335 GetESDsData(23)->Fill(towZNA_lg[i]);
336 GetESDsData(24)->Fill(towZPC_lg[i]);
337 GetESDsData(25)->Fill(towZPA_lg[i]);*/
340 sumQZNC += towZNC[i];
341 sumQZPC += towZPC[i];
342 sumQZNA += towZNA[i];
343 sumQZPA += towZPA[i];
345 /*sumQZNC_lg += towZNC_lg[i];
346 sumQZPC_lg += towZPC_lg[i];
347 sumQZNA_lg += towZNA_lg[i];
348 sumQZPA_lg += towZPA_lg[i];*/
351 GetESDsData(6)->Fill(sumQZNC);
352 GetESDsData(7)->Fill(sumQZNA);
353 GetESDsData(8)->Fill(sumQZPC);
354 GetESDsData(9)->Fill(sumQZPA);
356 /*GetESDsData(18)->Fill(sumQZNC_lg);
357 GetESDsData(19)->Fill(sumQZNA_lg);
358 GetESDsData(20)->Fill(sumQZPC_lg);
359 GetESDsData(21)->Fill(sumQZPA_lg);*/
362 //____________________________________________________________________________
363 void AliZDCQADataMakerRec::StartOfDetectorCycle()
365 //Detector specific actions at start of cycle
369 //____________________________________________________________________________
370 void AliZDCQADataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX_t task, TObjArray * list)
372 //Detector specific actions at end of cycle
373 // do the QA checking
374 AliQAChecker::Instance()->Run(AliQA::kZDC, task, list) ;