Changes in expert 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 <TProfile.h>
30 #include <Riostream.h>
31 // --- Standard library ---
32
33 // --- AliRoot header files ---
34 #include "AliLog.h"
35 #include "AliQAChecker.h"
36 #include "AliZDCReco.h"
37 #include "AliRawReader.h"
38 #include "AliZDCQADataMakerRec.h"
39 #include "AliZDCRawStream.h"
40 #include "AliZDCDigit.h"
41 #include "AliESDZDC.h"
42 #include "AliESDEvent.h"
43
44 ClassImp(AliZDCQADataMakerRec)
45            
46 //____________________________________________________________________________ 
47   AliZDCQADataMakerRec::AliZDCQADataMakerRec() : 
48   AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kZDC), "ZDC Quality Assurance Data Maker"), 
49   fDigit(0)
50 {
51   // ctor
52 }
53
54 //____________________________________________________________________________ 
55 AliZDCQADataMakerRec::AliZDCQADataMakerRec(const AliZDCQADataMakerRec& qadm) :
56   AliQADataMakerRec(),      
57   fDigit(0)
58 {
59   //copy ctor 
60   SetName((const char*)qadm.GetName()); 
61   SetTitle((const char*)qadm.GetTitle()); 
62 }
63
64 //__________________________________________________________________
65 AliZDCQADataMakerRec& AliZDCQADataMakerRec::operator = (const AliZDCQADataMakerRec& qadm )
66 {
67   // Equal operator.
68   this->~AliZDCQADataMakerRec();
69   new(this) AliZDCQADataMakerRec(qadm);
70   return *this;
71 }
72
73 //____________________________________________________________________________ 
74 void AliZDCQADataMakerRec::InitDigits()
75 {
76   // create Digits histograms in Digits subdir
77   //
78   const Bool_t expert   = kTRUE ; 
79   const Bool_t image    = kTRUE ; 
80   
81   // ------------------- HIGH GAIN CHAIN ---------------------------
82   TH1F * hDigZNCTot = new TH1F("hDigZNCTot", "Signal in ZNC;Amplitude [ADC counts];Counts", 100, 0., 6000.);
83   TH1F * hDigZNATot = new TH1F("hDigZNATot", "Signal in ZNA;Amplitude [ADC counts];Counts", 100, 0., 6000.);
84   TH1F * hDigZPCTot = new TH1F("hDigZPCTot", "Signal in ZPC;Amplitude [ADC counts];Counts", 100, 0., 6000.);
85   TH1F * hDigZPATot = new TH1F("hDigZPATot", "Signal in ZPA;Amplitude [ADC counts];Counts", 100, 0., 6000.);
86   Add2DigitsList(hDigZNCTot, 0, !expert, image);
87   Add2DigitsList(hDigZNATot, 1, !expert, image);
88   Add2DigitsList(hDigZPCTot, 2, !expert, image);
89   Add2DigitsList(hDigZPATot, 3, !expert, image);
90   //
91   TH1F * hDigSumQZNC = new TH1F("hDigSumQZNC", "Signal in 4 ZNC PMQ;Amplitude [ADC counts];Counts",100, 0., 4000.);
92   TH1F * hDigSumQZNA = new TH1F("hDigSumQZNA", "Signal in 4 ZNA PMQ;Amplitude [ADC counts];Counts",100, 0., 4000.);
93   TH1F * hDigSumQZPC = new TH1F("hDigSumQZPC", "Signal in 4 ZPC PMQ;Amplitude [ADC counts];Counts",100, 0., 4000.);
94   TH1F * hDigSumQZPA = new TH1F("hDigSumQZPA", "Signal in 4 ZPA PMQ;Amplitude [ADC counts];Counts",100, 0., 4000.);
95   Add2DigitsList(hDigSumQZNC, 4, expert, !image);
96   Add2DigitsList(hDigSumQZNA, 5, expert, !image);
97   Add2DigitsList(hDigSumQZPC, 6, expert, !image);
98   Add2DigitsList(hDigSumQZPA, 7, expert, !image);
99   //
100   TH1F * hDigPMCZNC = new TH1F("hDigPMCZNC", "Signal in ZNC PMC;Amplitude [ADC counts];Counts",100, 0., 4000.);
101   TH1F * hDigPMCZNA = new TH1F("hDigPMCZNA", "Signal in ZNA PMC;Amplitude [ADC counts];Counts",100, 0., 4000.);
102   TH1F * hDigPMCZPC = new TH1F("hDigPMCZPC", "Signal in ZPC PMC;Amplitude [ADC counts];Counts",100, 0., 4000.);
103   TH1F * hDigPMCZPA = new TH1F("hDigPMCZPA", "Signal in ZPA PMC;Amplitude [ADC counts];Counts",100, 0., 4000.);
104   Add2DigitsList(hDigPMCZNC, 8, expert, !image);
105   Add2DigitsList(hDigPMCZNA, 9, expert, !image);
106   Add2DigitsList(hDigPMCZPC, 10, expert, !image);
107   Add2DigitsList(hDigPMCZPA, 11, expert, !image);
108   // 
109   // ------------------- LOW GAIN CHAIN ---------------------------
110   TH1F * hDigZNCTotlg = new TH1F("hDigZNCTotlg", "Digit lg signal in ZNC", 100, 0., 6000.);
111   TH1F * hDigZNATotlg = new TH1F("hDigZNATotlg", "Digit lg signal in ZNA", 100, 0., 6000.);
112   TH1F * hDigZPCTotlg = new TH1F("hDigZPCTotlg", "Digit lg signal in ZPC", 100, 0., 6000.);
113   TH1F * hDigZPATotlg = new TH1F("hDigZPATotlg", "Digit lg signal in ZPA", 100, 0., 6000.);
114   Add2DigitsList(hDigZNCTotlg, 12, !expert, image);
115   Add2DigitsList(hDigZNATotlg, 13, !expert, image);
116   Add2DigitsList(hDigZPCTotlg, 14, !expert, image);
117   Add2DigitsList(hDigZPATotlg, 15, !expert, image);
118   //
119   TH1F * hDigSumQZNClg = new TH1F("hDigSumQZNClg", "Signal in 4 ZNC PMQlg",100, 0., 4000.);
120   TH1F * hDigSumQZNAlg = new TH1F("hDigSumQZNAlg", "Signal in 4 ZNA PMQlg",100, 0., 4000.);
121   TH1F * hDigSumQZPClg = new TH1F("hDigSumQZPClg", "Signal in 4 ZPC PMQlg",100, 0., 4000.);
122   TH1F * hDigSumQZPAlg = new TH1F("hDigSumQZPAlg", "Signal in 4 ZPA PMQlg",100, 0., 4000.);
123   Add2DigitsList(hDigSumQZNClg, 16, expert, !image);
124   Add2DigitsList(hDigSumQZNAlg, 17, expert, !image);
125   Add2DigitsList(hDigSumQZPClg, 18, expert, !image);
126   Add2DigitsList(hDigSumQZPAlg, 19, expert, !image);
127   //
128   TH1F * hDigPMCZNClg = new TH1F("hDigPMCZNClg", "Signal in ZNC PMClg",100, 0., 4000.);
129   TH1F * hDigPMCZNAlg = new TH1F("hDigPMCZNAlg", "Signal in ZNA PMClg",100, 0., 4000.);
130   TH1F * hDigPMCZPClg = new TH1F("hDigPMCZPClg", "Signal in ZPC PMClg",100, 0., 4000.);
131   TH1F * hDigPMCZPAlg = new TH1F("hDigPMCZPAlg", "Signal in ZPA PMClg",100, 0., 4000.);
132   Add2DigitsList(hDigPMCZNClg, 20, expert, !image);
133   Add2DigitsList(hDigPMCZNAlg, 21, expert, !image);
134   Add2DigitsList(hDigPMCZPClg, 22, expert, !image);
135   Add2DigitsList(hDigPMCZPAlg, 23, expert, !image);
136   
137 }
138
139 //____________________________________________________________________________
140 void AliZDCQADataMakerRec::InitRecPoints()
141 {
142   // create Digits histograms in Digits subdir
143   const Bool_t expert = kTRUE ; 
144   const Bool_t image  = kTRUE ; 
145   //
146   // ------------------- HIGH GAIN CHAIN ---------------------------
147   TH1F * hRecZNCTot = new TH1F("hRecZNCTot", "Rec signal in ZNC;Amplitude [ADC counts];Counts", 100, 0., 6000.);
148   TH1F * hRecZNATot = new TH1F("hRecZNATot", "Rec signal in ZNA;Amplitude [ADC counts];Counts", 100, 0., 6000.);
149   TH1F * hRecZPCTot = new TH1F("hRecZPCTot", "Rec signal in ZPC;Amplitude [ADC counts];Counts", 100, 0., 10000.);
150   TH1F * hRecZPATot = new TH1F("hRecZPATot", "Rec signal in ZPA;Amplitude [ADC counts];Counts", 100, 0., 10000.);
151   Add2RecPointsList(hRecZNCTot, 0, expert, !image);
152   Add2RecPointsList(hRecZNATot, 1, expert, !image);
153   Add2RecPointsList(hRecZPCTot, 2, expert, !image);
154   Add2RecPointsList(hRecZPATot, 3, expert, !image);
155   //
156   TH1F * hRecSumQZNC = new TH1F("hRecSumQZNC", "Rec summed 4 ZNC quadrants;Amplitude [ADC counts];Counts",100, 0., 4000.);
157   TH1F * hRecSumQZNA = new TH1F("hRecSumQZNA", "Rec summed 4 ZNA quadrants;Amplitude [ADC counts];Counts",100, 0., 4000.);
158   TH1F * hRecSumQZPC = new TH1F("hRecSumQZPC", "Rec summed 4 ZPC quadrants;Amplitude [ADC counts];Counts",100, 0., 4000.);
159   TH1F * hRecSumQZPA = new TH1F("hRecSumQZPA", "Rec summed 4 ZPA quadrants;Amplitude [ADC counts];Counts",100, 0., 4000.);
160   Add2RecPointsList(hRecSumQZNC, 4, expert, !image);
161   Add2RecPointsList(hRecSumQZNA, 5, expert, !image);
162   Add2RecPointsList(hRecSumQZPC, 6, expert, !image);
163   Add2RecPointsList(hRecSumQZPA, 7, expert, !image);
164   //
165   TH1F * hRecPMCZNC = new TH1F("hRecPMCZNC", "Rec common ZNC PMT;Amplitude [ADC counts];Counts",100, 0., 4000.);
166   TH1F * hRecPMCZNA = new TH1F("hRecPMCZNA", "Rec common ZNA PMT;Amplitude [ADC counts];Counts",100, 0., 4000.);
167   TH1F * hRecPMCZPC = new TH1F("hRecPMCZPC", "Rec common ZPC PMT;Amplitude [ADC counts];Counts",100, 0., 4000.);
168   TH1F * hRecPMCZPA = new TH1F("hRecPMCZPA", "Rec common ZPA PMT;Amplitude [ADC counts];Counts",100, 0., 4000.);
169   Add2RecPointsList(hRecPMCZNC, 8 , !expert, image);
170   Add2RecPointsList(hRecPMCZNA, 9 , !expert, image);
171   Add2RecPointsList(hRecPMCZPC, 10, !expert, image);
172   Add2RecPointsList(hRecPMCZPA, 11, !expert, image); 
173 }
174
175
176 //____________________________________________________________________________
177 void AliZDCQADataMakerRec::InitRaws()
178 {
179   // create Digits histograms in Digits subdir
180   const Bool_t expert   = kTRUE ; 
181   const Bool_t image    = kTRUE ; 
182   //
183   // ------------------- HIGH GAIN CHAIN ---------------------------
184   TH1F * hRawZNCTot = new TH1F("hRawZNCTot", "Raw signal in ZNC;Amplitude [ADC counts];Counts", 100, 0., 6000.);
185   TH1F * hRawZNATot = new TH1F("hRawZNATot", "Raw signal in ZNA;Amplitude [ADC counts];Counts", 100, 0., 6000.);
186   TH1F * hRawZPCTot = new TH1F("hRawZPCTot", "Raw signal in ZPC;Amplitude [ADC counts];Counts", 100, 0., 10000.);
187   TH1F * hRawZPATot = new TH1F("hRawZPATot", "Raw signal in ZPA;Amplitude [ADC counts];Counts", 100, 0., 10000.);
188   Add2RawsList(hRawZNCTot, 0, expert, !image);
189   Add2RawsList(hRawZNATot, 1, expert, !image);
190   Add2RawsList(hRawZPCTot, 2, expert, !image);
191   Add2RawsList(hRawZPATot, 3, expert, !image);
192   //
193   TH1F * hRawSumQZNC = new TH1F("hRawSumQZNC", "Raw summed 4 ZNC quadrants;Amplitude [ADC counts];Counts",100, 0., 4000.);
194   TH1F * hRawSumQZNA = new TH1F("hRawSumQZNA", "Raw summed 4 ZNA quadrants;Amplitude [ADC counts];Counts",100, 0., 4000.);
195   TH1F * hRawSumQZPC = new TH1F("hRawSumQZPC", "Raw summed 4 ZPC quadrants;Amplitude [ADC counts];Counts",100, 0., 4000.);
196   TH1F * hRawSumQZPA = new TH1F("hRawSumQZPA", "Raw summed 4 ZPA quadrants;Amplitude [ADC counts];Counts",100, 0., 4000.);
197   Add2RawsList(hRawSumQZNC, 4, expert, !image);
198   Add2RawsList(hRawSumQZNA, 5, expert, !image);
199   Add2RawsList(hRawSumQZPC, 6, expert, !image);
200   Add2RawsList(hRawSumQZPA, 7, expert, !image);
201   //
202   TH1F * hRawPMCZNC = new TH1F("hRawPMCZNC", "Raw common ZNC PMT;Amplitude [ADC counts];Counts",100, 0., 4000.);
203   TH1F * hRawPMCZNA = new TH1F("hRawPMCZNA", "Raw common ZNA PMT;Amplitude [ADC counts];Counts",100, 0., 4000.);
204   TH1F * hRawPMCZPC = new TH1F("hRawPMCZPC", "Raw common ZPC PMT;Amplitude [ADC counts];Counts",100, 0., 4000.);
205   TH1F * hRawPMCZPA = new TH1F("hRawPMCZPA", "Raw common ZPA PMT;Amplitude [ADC counts];Counts",100, 0., 4000.);
206   Add2RawsList(hRawPMCZNC, 8 , !expert, image);
207   Add2RawsList(hRawPMCZNA, 9 , !expert, image);
208   Add2RawsList(hRawPMCZPC, 10, !expert, image);
209   Add2RawsList(hRawPMCZPA, 11, !expert, image);
210   // 
211   // ------------------- LOW GAIN CHAIN ---------------------------
212   TH1F * hRawZNCTotlg = new TH1F("hRawZNCTotlg", "Rawit lg signal in ZNC", 100, 0., 6000.);
213   TH1F * hRawZNATotlg = new TH1F("hRawZNATotlg", "Rawit lg signal in ZNA", 100, 0., 6000.);
214   TH1F * hRawZPCTotlg = new TH1F("hRawZPCTotlg", "Rawit lg signal in ZPC", 100, 0., 10000.);
215   TH1F * hRawZPATotlg = new TH1F("hRawZPATotlg", "Rawit lg signal in ZPA", 100, 0., 10000.);
216   Add2RawsList(hRawZNCTotlg, 12, expert, !image);
217   Add2RawsList(hRawZNATotlg, 13, expert, !image);
218   Add2RawsList(hRawZPCTotlg, 14, expert, !image);
219   Add2RawsList(hRawZPATotlg, 15, expert, !image);
220   //
221   TH1F * hRawSumQZNClg = new TH1F("hRawSumQZNClg", "Raw summed 4 lg ZNC quadrants",100, 0., 4000.);
222   TH1F * hRawSumQZNAlg = new TH1F("hRawSumQZNAlg", "Raw summed 4 lg ZNA quadrants",100, 0., 4000.);
223   TH1F * hRawSumQZPClg = new TH1F("hRawSumQZPClg", "Raw summed 4 lg ZPC quadrants",100, 0., 4000.);
224   TH1F * hRawSumQZPAlg = new TH1F("hRawSumQZPAlg", "Raw summed 4 lg ZPA quadrants",100, 0., 4000.);
225   Add2RawsList(hRawSumQZNClg, 16, expert, !image);
226   Add2RawsList(hRawSumQZNAlg, 17, expert, !image);
227   Add2RawsList(hRawSumQZPClg, 18, expert, !image);
228   Add2RawsList(hRawSumQZPAlg, 19, expert, !image);
229   //
230   TH1F * hRawPMCZNClg = new TH1F("hRawPMCZNClg", "Raw common lg ZNC PMT",100, 0., 4000.);
231   TH1F * hRawPMCZNAlg = new TH1F("hRawPMCZNAlg", "Raw common lg ZNA PMT",100, 0., 4000.);
232   TH1F * hRawPMCZPClg = new TH1F("hRawPMCZPClg", "Raw common lg ZPC PMT",100, 0., 4000.);
233   TH1F * hRawPMCZPAlg = new TH1F("hRawPMCZPAlg", "Raw common lg ZPA PMT",100, 0., 4000.);
234   Add2RawsList(hRawPMCZNClg, 20, expert, !image);
235   Add2RawsList(hRawPMCZNAlg, 21, expert, !image);
236   Add2RawsList(hRawPMCZPClg, 22, expert, !image);
237   Add2RawsList(hRawPMCZPAlg, 23, expert, !image);
238 }
239
240 //____________________________________________________________________________
241 void AliZDCQADataMakerRec::InitESDs()
242 {
243   //Booking ESDs histograms
244   //
245   const Bool_t expert   = kTRUE ; 
246   const Bool_t image    = kTRUE ; 
247   
248   TH2F * hZNC  = new TH2F("hZNC", "Centroid in ZNC", 100, -5.,5.,100,-5.,5.);
249   TH2F * hZNA  = new TH2F("hZNA", "Centroid in ZNA", 100, -5.,5.,100,-5.,5.);
250   Add2ESDsList(hZNC, 0, !expert, image);
251   Add2ESDsList(hZNA, 1, !expert, image);
252   //
253   // ------------------- HIGH GAIN CHAIN ---------------------------
254   TH1F * hESDZNCTot = new TH1F("hESDZNCTot", "Energy in ZNC", 100, 0., 6000.);
255   TH1F * hESDZNATot = new TH1F("hESDZNATot", "Energy in ZNA", 100, 0., 6000.);
256   TH1F * hESDZPCTot = new TH1F("hESDZPCTot", "Energy in ZPC", 100, 0., 10000.);
257   TH1F * hESDZPATot = new TH1F("hESDZPATot", "Energy in ZPA", 100, 0., 10000.);
258   Add2ESDsList(hESDZNCTot, 2, expert, !image);
259   Add2ESDsList(hESDZNATot, 3, expert, !image);
260   Add2ESDsList(hESDZPCTot, 4, expert, !image);
261   Add2ESDsList(hESDZPATot, 5, expert, !image);
262   //
263   TH1F * hESDSumQZNC = new TH1F("hESDSumQZNC", "Sum of 4 ZNC energy",100, 0., 4000.);
264   TH1F * hESDSumQZNA = new TH1F("hESDSumQZNA", "Sum of 4 ZNA energy",100, 0., 4000.);
265   TH1F * hESDSumQZPC = new TH1F("hESDSumQZPC", "Sum of 4 ZPC energy",100, 0., 4000.);
266   TH1F * hESDSumQZPA = new TH1F("hESDSumQZPA", "Sum of 4 ZPA energy",100, 0., 4000.);
267   Add2ESDsList(hESDSumQZNC, 6, expert, !image);
268   Add2ESDsList(hESDSumQZNA, 7, expert, !image);
269   Add2ESDsList(hESDSumQZPC, 8, expert, !image);
270   Add2ESDsList(hESDSumQZPA, 9, expert, !image);
271   //
272   TH1F * hESDPMCZNC = new TH1F("hESDPMCZNC", "Energy in common ZNC PMT",100, 0., 4000.);
273   TH1F * hESDPMCZNA = new TH1F("hESDPMCZNA", "Energy in common ZNA PMT",100, 0., 4000.);
274   TH1F * hESDPMCZPC = new TH1F("hESDPMCZPC", "Energy in common ZPC PMT",100, 0., 4000.);
275   TH1F * hESDPMCZPA = new TH1F("hESDPMCZPA", "Energy in common ZPA PMT",100, 0., 4000.);
276   Add2ESDsList(hESDPMCZNC, 10, !expert, image);
277   Add2ESDsList(hESDPMCZNA, 11, !expert, image);
278   Add2ESDsList(hESDPMCZPC, 12, !expert, image);
279   Add2ESDsList(hESDPMCZPA, 13, !expert, image);
280   // 
281   // ------------------- LOW GAIN CHAIN ---------------------------
282   TH1F * hESDZNCTotlg = new TH1F("hESDZNCTotlg", "ESD lg signal in ZNC", 100, 0., 6000.);
283   TH1F * hESDZNATotlg = new TH1F("hESDZNATotlg", "ESD lg signal in ZNA", 100, 0., 6000.);
284   TH1F * hESDZPCTotlg = new TH1F("hESDZPCTotlg", "ESD lg signal in ZPC", 100, 0., 10000.);
285   TH1F * hESDZPATotlg = new TH1F("hESDZPATotlg", "ESD lg signal in ZPA", 100, 0., 10000.);
286   Add2ESDsList(hESDZNCTotlg, expert, !image);
287   Add2ESDsList(hESDZNATotlg, expert, !image);
288   Add2ESDsList(hESDZPCTotlg, expert, !image);
289   Add2ESDsList(hESDZPATotlg, expert, !image);
290   //
291   TH1F * hESDSumQZNClg = new TH1F("hESDSumQZNClg", "Sum of 4 lg ZNC sectors",100, 0., 4000.);
292   TH1F * hESDSumQZNAlg = new TH1F("hESDSumQZNAlg", "Sum of 4 lg ZNA sectors",100, 0., 4000.);
293   TH1F * hESDSumQZPClg = new TH1F("hESDSumQZPClg", "Sum of 4 lg ZPC sectors",100, 0., 4000.);
294   TH1F * hESDSumQZPAlg = new TH1F("hESDSumQZPAlg", "Sum of 4 lg ZPA sectors",100, 0., 4000.);
295   Add2ESDsList(hESDSumQZNClg, 18, expert, !image);
296   Add2ESDsList(hESDSumQZNAlg, 19, expert, !image);
297   Add2ESDsList(hESDSumQZPClg, 20, expert, !image);
298   Add2ESDsList(hESDSumQZPAlg, 21, expert, !image);
299   //
300   TH1F * hESDPMCZNClg = new TH1F("hESDPMCZNClg", "Signal in common ZNC lg PMT",100, 0., 4000.);
301   TH1F * hESDPMCZNAlg = new TH1F("hESDPMCZNAlg", "Signal in common ZNA lg PMT",100, 0., 4000.);
302   TH1F * hESDPMCZPClg = new TH1F("hESDPMCZPClg", "Signal in common ZPC lg PMT",100, 0., 4000.);
303   TH1F * hESDPMCZPAlg = new TH1F("hESDPMCZPAlg", "Signal in common ZPA lg PMT",100, 0., 4000.);
304   Add2ESDsList(hESDPMCZNClg, 22, expert, !image);
305   Add2ESDsList(hESDPMCZNAlg, 23, expert, !image);
306   Add2ESDsList(hESDPMCZPClg, 24, expert, !image);
307   Add2ESDsList(hESDPMCZPAlg, 25, expert, !image);
308 }
309
310 //___________________________________________________________________________
311 void AliZDCQADataMakerRec::MakeDigits(TTree *digitTree)
312 {
313   TBranch * branch = digitTree->GetBranch("ZDC");
314   if(!branch){
315     AliError("ZDC branch in Digit Tree not found"); 
316     return;
317   } 
318   
319   // Check id histograms already created for this Event Specie
320   if(!GetDigitsData(0)) InitDigits() ;
321   
322   branch->SetAddress(&fDigit);
323   
324   Int_t ndig = digitTree->GetEntries();
325   
326   Float_t adcSum_ZNC=0., adcSum_ZNA=0., adcSum_ZPC=0., adcSum_ZPA=0.;
327   Float_t adcSumQ_ZNC=0., adcSumQ_ZNA=0., adcSumQ_ZPC=0., adcSumQ_ZPA=0.;
328   Float_t adcSum_ZNC_lg=0., adcSum_ZNA_lg=0., adcSum_ZPC_lg=0., adcSum_ZPA_lg=0.;
329   Float_t adcSumQ_ZNC_lg=0., adcSumQ_ZNA_lg=0., adcSumQ_ZPC_lg=0., adcSumQ_ZPA_lg=0.;
330   //
331   for(Int_t i = 0; i < ndig; i++){
332     digitTree->GetEntry(i);
333     if(fDigit->GetSector(0)==1){
334       adcSum_ZNC += fDigit->GetADCValue(0);
335       adcSum_ZNC_lg += fDigit->GetADCValue(1);
336       //
337       if(fDigit->GetSector(1)!=0){
338               adcSumQ_ZNC += fDigit->GetADCValue(0);
339               adcSumQ_ZNC_lg+= fDigit->GetADCValue(1);
340       }
341       else{
342               GetDigitsData(8)->Fill(fDigit->GetADCValue(0));
343               GetDigitsData(20)->Fill(fDigit->GetADCValue(1));
344       }
345     }
346     else if(fDigit->GetSector(0)==2){
347       adcSum_ZPC += fDigit->GetADCValue(0);
348       adcSum_ZPC_lg += fDigit->GetADCValue(1);
349       //
350       if(fDigit->GetSector(1)!=0){
351               adcSumQ_ZPC += fDigit->GetADCValue(0);
352               adcSumQ_ZPC_lg+= fDigit->GetADCValue(1);
353       }
354       else{
355               GetDigitsData(10)->Fill(fDigit->GetADCValue(0));
356               GetDigitsData(22)->Fill(fDigit->GetADCValue(1));
357       }
358     }
359     else if(fDigit->GetSector(0)==4){
360       adcSum_ZNA += fDigit->GetADCValue(0);
361       adcSum_ZNA_lg += fDigit->GetADCValue(1);
362       //
363       if(fDigit->GetSector(1)!=0){
364               adcSumQ_ZNA += fDigit->GetADCValue(0);
365               adcSumQ_ZNA_lg+= fDigit->GetADCValue(1);
366       }
367       else{
368               GetDigitsData(9)->Fill(fDigit->GetADCValue(0));
369               GetDigitsData(21)->Fill(fDigit->GetADCValue(1));
370       }
371     }
372     else if(fDigit->GetSector(0)==5){
373       adcSum_ZPA += fDigit->GetADCValue(0);
374       adcSum_ZPA_lg += fDigit->GetADCValue(1);
375       //
376       if(fDigit->GetSector(1)!=0){
377               adcSumQ_ZPA += fDigit->GetADCValue(0);
378               adcSumQ_ZPA_lg+= fDigit->GetADCValue(1);
379       }
380       else{
381               GetDigitsData(11)->Fill(fDigit->GetADCValue(0));
382               GetDigitsData(23)->Fill(fDigit->GetADCValue(1));
383       }
384     }
385   }
386   //
387   GetDigitsData(0)->Fill(adcSum_ZNC);
388   GetDigitsData(1)->Fill(adcSum_ZNA);
389   GetDigitsData(2)->Fill(adcSum_ZPC);
390   GetDigitsData(3)->Fill(adcSum_ZPA);
391   //
392   GetDigitsData(4)->Fill(adcSumQ_ZNC);
393   GetDigitsData(5)->Fill(adcSumQ_ZNA);
394   GetDigitsData(6)->Fill(adcSumQ_ZPC);
395   GetDigitsData(7)->Fill(adcSumQ_ZPA);
396   //
397   GetDigitsData(12)->Fill(adcSum_ZNC_lg);
398   GetDigitsData(13)->Fill(adcSum_ZNA_lg);
399   GetDigitsData(14)->Fill(adcSum_ZPC_lg);
400   GetDigitsData(15)->Fill(adcSum_ZPA_lg);
401   //
402   GetDigitsData(16)->Fill(adcSumQ_ZNC_lg);
403   GetDigitsData(17)->Fill(adcSumQ_ZNA_lg);
404   GetDigitsData(18)->Fill(adcSumQ_ZPC_lg);
405   GetDigitsData(19)->Fill(adcSumQ_ZPA_lg);
406 }
407
408 //____________________________________________________________________________
409 void AliZDCQADataMakerRec::MakeRecPoints(TTree * clustersTree)
410 {
411   // Filling QA histos from RecPoints
412
413   TBranch *branch = clustersTree->GetBranch("ZDC");
414   if(!branch){ 
415     AliError("Can't get the ZDC branch for rec points!");
416     return;
417   }
418   
419   if(!GetRecPointsData(0)) InitRecPoints() ;
420
421   Float_t sum_ZNC=0., sum_ZNA=0., sum_ZPC=0., sum_ZPA=0.;
422   Float_t sumQ_ZNC=0., sumQ_ZNA=0., sumQ_ZPC=0., sumQ_ZPA=0.;
423
424   AliZDCReco reco;
425   AliZDCReco* preco = &reco;
426   clustersTree->SetBranchAddress("ZDC", &preco);
427
428   clustersTree->GetEntry(0);
429   for(Int_t i=0; i<5; i++){
430     sum_ZNC += reco.GetZN1HREnTow(i);
431     sum_ZPC += reco.GetZN2HREnTow(i);
432     sum_ZNA += reco.GetZP1HREnTow(i);
433     sum_ZPA += reco.GetZP2HREnTow(i);
434     if(i==0){
435       GetRecPointsData(8)->Fill(reco.GetZN1HREnTow(i));
436       GetRecPointsData(9)->Fill(reco.GetZN2HREnTow(i));
437       GetRecPointsData(10)->Fill(reco.GetZP1HREnTow(i));
438       GetRecPointsData(11)->Fill(reco.GetZP2HREnTow(i));
439     }
440     else{
441       sumQ_ZNC += reco.GetZN1HREnTow(i);
442       sumQ_ZPC += reco.GetZN2HREnTow(i);
443       sumQ_ZNA += reco.GetZP1HREnTow(i);
444       sumQ_ZPA += reco.GetZP2HREnTow(i);
445     }
446   }
447   
448   GetRecPointsData(0)->Fill(sum_ZNC);
449   GetRecPointsData(1)->Fill(sum_ZPC);
450   GetRecPointsData(2)->Fill(sum_ZNA);
451   GetRecPointsData(3)->Fill(sum_ZPA);
452   //
453   GetRecPointsData(4)->Fill(sumQ_ZNC);
454   GetRecPointsData(5)->Fill(sumQ_ZPC);
455   GetRecPointsData(6)->Fill(sumQ_ZNA);
456   GetRecPointsData(7)->Fill(sumQ_ZPA);
457   
458 }  
459
460 //____________________________________________________________________________
461 void AliZDCQADataMakerRec::MakeRaws(AliRawReader *rawReader)
462 {
463   // Filling Raws QA histos
464   //
465   // Check if histograms already created for this Event Specie
466   if(!GetRawsData(0)) InitRaws();
467
468   Float_t sum_ZNC=0., sum_ZNA=0., sum_ZPC=0., sum_ZPA=0.;
469   Float_t sumQ_ZNC=0., sumQ_ZNA=0., sumQ_ZPC=0., sumQ_ZPA=0.;
470   Float_t sum_ZNC_lg=0., sum_ZNA_lg=0., sum_ZPC_lg=0., sum_ZPA_lg=0.;
471   Float_t sumQ_ZNC_lg=0., sumQ_ZNA_lg=0., sumQ_ZPC_lg=0., sumQ_ZPA_lg=0.;
472   //
473   AliZDCRawStream stream(rawReader);
474   while(stream.Next()){
475     if(stream.IsADCDataWord() && 
476      (stream.GetADCModule()==0 || stream.GetADCModule()==1)){
477        if(stream.GetSector(0)==1){
478          if(stream.GetADCGain()==0){
479            sum_ZNC += stream.GetADCValue();
480            if(stream.GetSector(1)!=0) sumQ_ZNC += stream.GetADCValue();
481            else GetRawsData(8)->Fill(stream.GetADCValue());
482          }
483          else{
484            sum_ZNC_lg += stream.GetADCValue();
485            if(stream.GetSector(1)!=0) sumQ_ZNC_lg += stream.GetADCValue();
486            else GetRawsData(20)->Fill(stream.GetADCValue());
487          }
488        }
489        else if(stream.GetSector(0)==2){
490          if(stream.GetADCGain()==0){
491            sum_ZPC += stream.GetADCValue();
492            if(stream.GetSector(1)!=0) sumQ_ZPC += stream.GetADCValue();
493            else GetRawsData(10)->Fill(stream.GetADCValue());
494          }
495          else{
496            sum_ZPC_lg += stream.GetADCValue();
497            if(stream.GetSector(1)!=0) sumQ_ZPC_lg += stream.GetADCValue();
498            else GetRawsData(22)->Fill(stream.GetADCValue());
499          }
500        }
501        else if(stream.GetSector(0)==4){
502          if(stream.GetADCGain()==0){
503            sum_ZNA += stream.GetADCValue();
504            if(stream.GetSector(1)!=0) sumQ_ZNA += stream.GetADCValue();
505            else GetRawsData(9)->Fill(stream.GetADCValue());
506          }
507          else{
508            sum_ZNA_lg += stream.GetADCValue();
509            if(stream.GetSector(1)!=0) sumQ_ZNA_lg += stream.GetADCValue();
510            else GetRawsData(21)->Fill(stream.GetADCValue());
511          }
512        }
513        else if(stream.GetSector(0)==5){
514          if(stream.GetADCGain()==0){
515            sum_ZPA += stream.GetADCValue();
516            if(stream.GetSector(1)!=0) sumQ_ZPA += stream.GetADCValue();
517            else GetRawsData(11)->Fill(stream.GetADCValue());
518          }
519          else{
520            sum_ZPA_lg += stream.GetADCValue();
521            if(stream.GetSector(1)!=0) sumQ_ZPA_lg += stream.GetADCValue();
522            else GetRawsData(23)->Fill(stream.GetADCValue());
523          }
524        }
525     }
526   }
527   //
528   GetRawsData(0)->Fill(sum_ZNC);
529   GetRawsData(1)->Fill(sum_ZNA);
530   GetRawsData(2)->Fill(sum_ZPC);
531   GetRawsData(3)->Fill(sum_ZPA);
532   //
533   GetRawsData(4)->Fill(sumQ_ZNC);
534   GetRawsData(5)->Fill(sumQ_ZNA);
535   GetRawsData(6)->Fill(sumQ_ZPC);
536   GetRawsData(7)->Fill(sumQ_ZPA);
537   //
538   GetRawsData(12)->Fill(sum_ZNC_lg);
539   GetRawsData(13)->Fill(sum_ZNA_lg);
540   GetRawsData(14)->Fill(sum_ZPC_lg);
541   GetRawsData(15)->Fill(sum_ZPA_lg);
542   //
543   GetRawsData(16)->Fill(sumQ_ZNC_lg);
544   GetRawsData(17)->Fill(sumQ_ZNA_lg);
545   GetRawsData(18)->Fill(sumQ_ZPC_lg);
546   GetRawsData(19)->Fill(sumQ_ZPA_lg);
547   //
548 //   stream.Delete();
549 }
550
551 //____________________________________________________________________________
552 void AliZDCQADataMakerRec::MakeESDs(AliESDEvent * esd)
553 {
554   // make QA data from ESDs
555   //
556   
557   // Check id histograms already created for this Event Specie
558   if ( ! GetESDsData(0) )
559     InitESDs() ;
560
561   AliESDZDC * zdcESD =  esd->GetESDZDC();
562   //
563   Double32_t * centr_ZNC = zdcESD->GetZNCCentroid();
564   GetESDsData(0)->Fill(centr_ZNC[0], centr_ZNC[1]);
565
566   Double32_t * centr_ZNA = zdcESD->GetZNACentroid();
567   GetESDsData(1)->Fill(centr_ZNA[0], centr_ZNA[1]);
568
569   //
570   GetESDsData(2)->Fill(esd->GetZDCN1Energy());
571   GetESDsData(3)->Fill(esd->GetZDCN2Energy());
572   GetESDsData(4)->Fill(esd->GetZDCP1Energy());
573   GetESDsData(5)->Fill(esd->GetZDCP2Energy());
574   //
575   Double_t sumQZNC=0., sumQZPC=0., sumQZNA=0., sumQZPA=0.;
576   Double_t sumQZNC_lg=0., sumQZPC_lg=0., sumQZNA_lg=0., sumQZPA_lg=0.;
577   //
578   const Double_t *towZNC, *towZPC, *towZNA, *towZPA;
579   const Double_t *towZNC_lg, *towZPC_lg, *towZNA_lg, *towZPA_lg;
580   //
581   towZNC = zdcESD->GetZN1TowerEnergy();
582   towZPC = zdcESD->GetZP1TowerEnergy();
583   towZNA = zdcESD->GetZN2TowerEnergy();
584   towZPA = zdcESD->GetZP2TowerEnergy();
585   //
586   towZNC_lg = zdcESD->GetZN1TowerEnergyLR();
587   towZPC_lg = zdcESD->GetZP1TowerEnergyLR();
588   towZNA_lg = zdcESD->GetZN2TowerEnergyLR();
589   towZPA_lg = zdcESD->GetZP2TowerEnergyLR();
590   //
591   for(Int_t i=0; i<5; i++){
592      if(i==0){
593        GetESDsData(10)->Fill(towZNC[i]);
594        GetESDsData(11)->Fill(towZNA[i]);
595        GetESDsData(12)->Fill(towZPC[i]);
596        GetESDsData(13)->Fill(towZPA[i]);
597        //
598        GetESDsData(22)->Fill(towZNC_lg[i]);
599        GetESDsData(23)->Fill(towZNA_lg[i]);
600        GetESDsData(24)->Fill(towZPC_lg[i]);
601        GetESDsData(25)->Fill(towZPA_lg[i]);
602      }
603      else{
604        sumQZNC += towZNC[i];
605        sumQZPC += towZPC[i];
606        sumQZNA += towZNA[i];
607        sumQZPA += towZPA[i];
608        //
609        sumQZNC_lg += towZNC_lg[i];
610        sumQZPC_lg += towZPC_lg[i];
611        sumQZNA_lg += towZNA_lg[i];
612        sumQZPA_lg += towZPA_lg[i];
613      }
614   }
615   GetESDsData(6)->Fill(sumQZNC);
616   GetESDsData(7)->Fill(sumQZNA);
617   GetESDsData(8)->Fill(sumQZPC);
618   GetESDsData(9)->Fill(sumQZPA);
619   //
620   GetESDsData(18)->Fill(sumQZNC_lg);
621   GetESDsData(19)->Fill(sumQZNA_lg);
622   GetESDsData(20)->Fill(sumQZPC_lg);
623   GetESDsData(21)->Fill(sumQZPA_lg);
624 }
625
626 //____________________________________________________________________________
627 void AliZDCQADataMakerRec::StartOfDetectorCycle()
628 {
629   //Detector specific actions at start of cycle
630   
631 }
632
633 //____________________________________________________________________________ 
634 void AliZDCQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
635 {
636   //Detector specific actions at end of cycle
637   // do the QA checking
638   AliQAChecker::Instance()->Run(AliQAv1::kZDC, task, list) ;  
639 }
640