]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSQADataMakerRec.cxx
Added QA for digits during reconstruction (Yves)
[u/mrichter/AliRoot.git] / PHOS / AliPHOSQADataMakerRec.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 /* $Id$ */
18
19 /*
20   Produces the data needed to calculate the quality assurance. 
21   All data must be mergeable objects.
22   Y. Schutz CERN July 2007
23 */
24
25 // --- ROOT system ---
26 #include <TClonesArray.h>
27 #include <TFile.h> 
28 #include <TH1F.h> 
29 #include <TH1I.h> 
30 #include <TH2F.h> 
31 #include <TParameter.h>
32
33 // --- Standard library ---
34
35 // --- AliRoot header files ---
36 #include "AliESDCaloCluster.h"
37 #include "AliESDEvent.h"
38 #include "AliLog.h"
39 #include "AliPHOSQADataMakerRec.h"
40 #include "AliQAChecker.h"
41 #include "AliPHOSDigit.h" 
42 #include "AliPHOSCpvRecPoint.h" 
43 #include "AliPHOSEmcRecPoint.h" 
44 #include "AliPHOSRecParticle.h" 
45 #include "AliPHOSTrackSegment.h" 
46 #include "AliPHOSRawDecoder.h"
47 #include "AliPHOSRawDecoderv1.h"
48 #include "AliPHOSRawDecoderv2.h"
49 #include "AliPHOSReconstructor.h"
50
51 ClassImp(AliPHOSQADataMakerRec)
52            
53 //____________________________________________________________________________ 
54   AliPHOSQADataMakerRec::AliPHOSQADataMakerRec() : 
55   AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kPHOS), "PHOS Quality Assurance Data Maker")
56 {
57   // ctor
58 }
59
60 //____________________________________________________________________________ 
61 AliPHOSQADataMakerRec::AliPHOSQADataMakerRec(const AliPHOSQADataMakerRec& qadm) :
62   AliQADataMakerRec()
63 {
64   //copy ctor 
65   SetName((const char*)qadm.GetName()) ; 
66   SetTitle((const char*)qadm.GetTitle()); 
67 }
68
69 //__________________________________________________________________
70 AliPHOSQADataMakerRec& AliPHOSQADataMakerRec::operator = (const AliPHOSQADataMakerRec& qadm )
71 {
72   // Equal operator.
73   this->~AliPHOSQADataMakerRec();
74   new(this) AliPHOSQADataMakerRec(qadm);
75   return *this;
76 }
77  
78 //____________________________________________________________________________ 
79 void AliPHOSQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
80 {
81   //Detector specific actions at end of cycle
82   if(GetRawsData(kHGqualMod1) && GetRawsData(kHGmod1))
83     GetRawsData(kHGqualMod1)->Divide( GetRawsData(kHGmod1) ) ;
84   if(GetRawsData(kHGqualMod2) && GetRawsData(kHGmod2))
85     GetRawsData(kHGqualMod2)->Divide( GetRawsData(kHGmod2) ) ;
86   if(GetRawsData(kHGqualMod3) && GetRawsData(kHGmod3))
87     GetRawsData(kHGqualMod3)->Divide( GetRawsData(kHGmod3) ) ;
88   if(GetRawsData(kHGqualMod4) && GetRawsData(kHGmod4))
89     GetRawsData(kHGqualMod4)->Divide( GetRawsData(kHGmod4) ) ;
90   if(GetRawsData(kHGqualMod5) && GetRawsData(kHGmod5))
91     GetRawsData(kHGqualMod5)->Divide( GetRawsData(kHGmod5) ) ;
92   // do the QA checking
93   AliQAChecker::Instance()->Run(AliQAv1::kPHOS, task, list) ;  
94 }
95
96 //____________________________________________________________________________ 
97 void AliPHOSQADataMakerRec::InitESDs()
98 {
99   //Create histograms to controll ESD
100  
101   const Bool_t expert   = kTRUE ; 
102   const Bool_t image    = kTRUE ; 
103
104   TH1F * h1 = new TH1F("hESDPhosSpectrum",  "ESDs spectrum in PHOS"                ,  200, 0.,   20.) ;
105   h1->Sumw2() ;
106   Add2ESDsList(h1, kESDSpec, !expert, image)  ;
107
108   TH1I * h2 = new TH1I("hESDPhosMul",       "ESDs multiplicity distribution in PHOS", 100, 0,   100 ) ; 
109   h2->Sumw2() ;
110   Add2ESDsList(h2, kESDNtot, !expert, image) ;
111  
112   TH1F * h3 = new TH1F("hESDPhosEtot",      "ESDs total energy"                     , 2000, 0,  200.) ; 
113   h3->Sumw2() ;
114   Add2ESDsList(h3, kESDEtot, !expert, image) ;  //Expert histo
115  
116   TH1F * h4 = new TH1F("hESDpid",           "ESDs PID distribution in PHOS"         , 100, 0.,    1.) ;
117   h4->Sumw2() ;
118   Add2ESDsList(h4, kESDpid, !expert, image) ; //Expert histo
119         
120 }
121
122 //____________________________________________________________________________ 
123 void AliPHOSQADataMakerRec::InitDigits()
124 {
125   // create Digits histograms in Digits subdir
126   const Bool_t expert   = kTRUE ; 
127   const Bool_t image    = kTRUE ; 
128   TH1I * h0 = new TH1I("hPhosDigits",    "Digits amplitude distribution in PHOS",    500, 0, 1000) ; 
129   h0->Sumw2() ;
130   Add2DigitsList(h0, kDigits, !expert, image) ;
131   TH1I * h1 = new TH1I("hPhosDigitsMul", "Digits multiplicity distribution in PHOS", 2000, 0, 10000) ; 
132   h1->Sumw2() ;
133   Add2DigitsList(h1, kDigitsMul, !expert, image) ;
134 }
135
136 //____________________________________________________________________________ 
137 void AliPHOSQADataMakerRec::InitRecPoints()
138 {
139   // create Reconstructed Points histograms in RecPoints subdir
140   const Bool_t expert   = kTRUE ; 
141   const Bool_t image    = kTRUE ; 
142   
143   TH2I * h0 = new TH2I("hRpPHOSxyMod1","RecPoints Rows x Columns for PHOS module 1", 64, -72., 72., 56, -63., 63.) ;                             
144   Add2RecPointsList(h0,kRPmod1, expert, !image) ;
145   TH2I * h1 = new TH2I("hRpPHOSxyMod2","RecPoints Rows x Columns for PHOS module 2", 64, -72., 72., 56, -63., 63.) ;                             
146   Add2RecPointsList(h1,kRPmod2, expert, !image) ;
147   TH2I * h2 = new TH2I("hRpPHOSxyMod3","RecPoints Rows x Columns for PHOS module 3", 64, -72., 72., 56, -63., 63.) ;                             
148   Add2RecPointsList(h2,kRPmod3, expert, !image) ;
149   TH2I * h3 = new TH2I("hRpPHOSxyMod4","RecPoints Rows x Columns for PHOS module 4", 64, -72., 72., 56, -63., 63.) ;                             
150   Add2RecPointsList(h3,kRPmod4, expert, !image) ;
151   TH2I * h4 = new TH2I("hRpPHOSxyMod5","RecPoints Rows x Columns for PHOS module 5", 64, -72., 72., 56, -63., 63.) ;                             
152   Add2RecPointsList(h4,kRPmod5, expert, !image) ;
153  
154   TH1F * h5 = new TH1F("hEmcPhosRecPointsSpectrum",  "EMC RecPoints spectrum in PHOS",   2000, 0., 20.) ; 
155   h5->Sumw2() ;
156   Add2RecPointsList(h5, kRPSpec, !expert, image)  ;
157
158   TH1I * h6 = new TH1I("hEmcPhosRecPointsMul", "EMCA RecPoints multiplicity distribution in PHOS", 100, 0,  100) ; 
159   h6->Sumw2() ;
160   Add2RecPointsList(h6, kRPNtot, !expert, image) ;
161
162   TH1I * h7 = new TH1I("hEmcPhosRecPointsEtot", "EMC RecPoints Etot", 200, 0,  200.) ; 
163   h7->Sumw2() ;
164   Add2RecPointsList(h7, kRPEtot, !expert, image) ;
165
166   TH1I * h8 = new TH1I("hCpvPhosRecPointsMul", "CPV RecPoints multiplicity distribution in PHOS", 100, 0,  100) ; 
167   h8->Sumw2() ;
168   Add2RecPointsList(h8, kRPNcpv, !expert, image) ;
169 }
170
171 //____________________________________________________________________________ 
172 void AliPHOSQADataMakerRec::InitRaws()
173 {
174   // create Raws histograms in Raws subdir
175   const Bool_t expert   = kTRUE ; 
176   const Bool_t saveCorr = kTRUE ; 
177   const Bool_t image    = kTRUE ; 
178
179   TH2I * h0 = new TH2I("hHighPHOSxyMod1","High Gain Rows x Columns for PHOS module 1", 64, 0, 64, 56, 0, 56) ;
180   Add2RawsList(h0,kHGmod1, expert, !image, !saveCorr) ;
181   TH2I * h1 = new TH2I("hHighPHOSxyMod2","High Gain Rows x Columns for PHOS module 2", 64, 0, 64, 56, 0, 56) ;
182   Add2RawsList(h1,kHGmod2, expert, !image, !saveCorr) ;
183   TH2I * h2 = new TH2I("hHighPHOSxyMod3","High Gain Rows x Columns for PHOS module 3", 64, 0, 64, 56, 0, 56) ;
184   Add2RawsList(h2,kHGmod3, expert, !image, !saveCorr) ;
185   TH2I * h3 = new TH2I("hHighPHOSxyMod4","High Gain Rows x Columns for PHOS module 4", 64, 0, 64, 56, 0, 56) ;
186   Add2RawsList(h3,kHGmod4, expert, !image, !saveCorr) ;
187   TH2I * h4 = new TH2I("hHighPHOSxyMod5","High Gain Rows x Columns for PHOS module 5", 64, 0, 64, 56, 0, 56) ;
188   Add2RawsList(h4,kHGmod5, expert, !image, !saveCorr) ;
189   TH2I * h5 = new TH2I("hLowPHOSxyMod1","Low Gain Rows x Columns for PHOS module 1", 64, 0, 64, 56, 0, 56) ;
190   Add2RawsList(h5,kLGmod1, expert, !image, !saveCorr) ;
191   TH2I * h6 = new TH2I("hLowPHOSxyMod2","Low Gain Rows x Columns for PHOS module 2", 64, 0, 64, 56, 0, 56) ;
192   Add2RawsList(h6,kLGmod2, expert, !image, !saveCorr) ;
193   TH2I * h7 = new TH2I("hLowPHOSxyMod3","Low Gain Rows x Columns for PHOS module 3", 64, 0, 64, 56, 0, 56) ;
194   Add2RawsList(h7,kLGmod3, expert, !image, !saveCorr) ;
195   TH2I * h8 = new TH2I("hLowPHOSxyMod4","Low Gain Rows x Columns for PHOS module 4", 64, 0, 64, 56, 0, 56) ;
196   Add2RawsList(h8,kLGmod4, expert, !image, !saveCorr) ;
197   TH2I * h9 = new TH2I("hLowPHOSxyMod5","Low Gain Rows x Columns for PHOS module 5", 64, 0, 64, 56, 0, 56) ;
198   Add2RawsList(h9,kLGmod5, expert, !image, !saveCorr) ;
199
200   TH1I * h10 = new TH1I("hLowPhosModules",    "Low Gain Hits in EMCA PHOS modules",       6, 0, 6) ;
201   h10->Sumw2() ;
202   Add2RawsList(h10, kNmodLG, !expert, image, !saveCorr) ;
203   TH1I * h11 = new TH1I("hHighPhosModules",   "High Gain Hits in EMCA PHOS modules",       6, 0, 6) ;
204   h11->Sumw2() ;
205   Add2RawsList(h11, kNmodHG, !expert, image, !saveCorr) ;
206
207   TH1F * h12 = new TH1F("hLowPhosRawtime", "Low Gain Time of raw hits in PHOS", 500, -50., 200.) ;
208   h12->Sumw2() ;
209   Add2RawsList(h12, kLGtime, !expert, image, !saveCorr) ;
210   TH1F * h13 = new TH1F("hHighPhosRawtime", "High Gain Time of raw hits in PHOS", 500, -50., 200.) ;
211   h13->Sumw2() ;
212   Add2RawsList(h13, kHGtime, !expert, image, !saveCorr) ;
213
214   TH1F * h14 = new TH1F("hLowPhosRawEnergy", "Low Gain Energy of raw hits in PHOS", 500, 0., 1000.) ;
215   h14->Sumw2() ;
216   Add2RawsList(h14, kSpecLG, !expert, image, !saveCorr) ;
217   TH1F * h15 = new TH1F("hHighPhosRawEnergy", "High Gain Energy of raw hits in PHOS",500,0., 1000.) ;
218   h15->Sumw2() ;
219   Add2RawsList(h15, kSpecHG, !expert, image, !saveCorr) ;
220
221   TH1F * h16 = new TH1F("hLowNtot", "Low Gain Total Number of raw hits in PHOS", 500, 0., 5000.) ;
222   h16->Sumw2() ;
223   Add2RawsList(h16, kNtotLG, !expert, saveCorr, image) ;
224   TH1F * h17 = new TH1F("hHighNtot", "High Gain Total Number of raw hits in PHOS",500,0., 5000.) ;
225   h17->Sumw2() ;
226   Add2RawsList(h17, kNtotHG, !expert, saveCorr, image) ;
227
228   TH1F * h18 = new TH1F("hLowEtot", "Low Gain Total Energy of raw hits in PHOS", 500, 0., 5000.) ;
229   h18->Sumw2() ;
230   Add2RawsList(h18, kEtotLG, !expert, saveCorr, image) ;
231   TH1F * h19 = new TH1F("hHighEtot", "High Gain Total Energy of raw hits in PHOS",500,0., 100000.) ;
232   h19->Sumw2() ;
233   Add2RawsList(h19, kEtotHG, !expert, saveCorr, image) ;
234
235   TH2F * h20 = new TH2F("hQualHGxyMod1","High Gain signal quality Rows x Columns module 1", 64, 0, 64, 56, 0, 56) ;
236   h20->SetOption("colz");
237   Add2RawsList(h20,kHGqualMod1, expert, !image, !saveCorr) ;
238   TH2F * h21 = new TH2F("hQualHGxyMod2","High Gain signal quality Rows x Columns module 2", 64, 0, 64, 56, 0, 56) ;
239   h21->SetOption("colz");
240   Add2RawsList(h21,kHGqualMod2, expert, !image, !saveCorr) ;
241   TH2F * h22 = new TH2F("hQualHGxyMod3","High Gain signal quality Rows x Columns module 3", 64, 0, 64, 56, 0, 56) ;
242   h22->SetOption("colz");
243   Add2RawsList(h22,kHGqualMod3, expert, !image, !saveCorr) ;
244   TH2F * h23 = new TH2F("hQualHGxyMod4","High Gain signal quality Rows x Columns module 4", 64, 0, 64, 56, 0, 56) ;
245   h23->SetOption("colz");
246   Add2RawsList(h23,kHGqualMod4, expert, !image, !saveCorr) ;
247   TH2F * h24 = new TH2F("hQualHGxyMod5","High Gain signal quality Rows x Columns module 5", 64, 0, 64, 56, 0, 56) ;
248   h24->SetOption("colz");
249   Add2RawsList(h24,kHGqualMod5, expert, !image, !saveCorr) ;
250   TH2F * h25 = new TH2F("hQualLGxyMod1","Low Gain signal quality Rows x Columns module 1", 64, 0, 64, 56, 0, 56) ;
251   h25->SetOption("colz");
252   Add2RawsList(h25,kLGqualMod1, expert, !image, !saveCorr) ;
253   TH2F * h26 = new TH2F("hQualLGxyMod2","Low Gain signal quality Rows x Columns module 2", 64, 0, 64, 56, 0, 56) ;
254   h26->SetOption("colz");
255   Add2RawsList(h26,kLGqualMod2, expert, !image, !saveCorr) ;
256   TH2F * h27 = new TH2F("hQualLGxyMod3","Low Gain signal quality Rows x Columns module 3", 64, 0, 64, 56, 0, 56) ;
257   h27->SetOption("colz");
258   Add2RawsList(h27,kLGqualMod3, expert, !image, !saveCorr) ;
259   TH2F * h28 = new TH2F("hQualLGxyMod4","Low Gain signal quality Rows x Columns module 4", 64, 0, 64, 56, 0, 56) ;
260   h28->SetOption("colz");
261   Add2RawsList(h28,kLGqualMod4, expert, !image, !saveCorr) ;
262   TH2F * h29 = new TH2F("hQualLGxyMod5","Low Gain signal quality Rows x Columns module 5", 64, 0, 64, 56, 0, 56) ;
263   h29->SetOption("colz");
264   Add2RawsList(h29,kLGqualMod5, expert, !image, !saveCorr) ;
265
266   TH1F * h30 = new TH1F("hLGpedRMS","Low Gain pedestal RMS",200,0.,20.) ;
267   h30->Sumw2() ;
268   Add2RawsList(h30,kLGpedRMS, !expert, image, !saveCorr) ;
269   TH1F * h31 = new TH1F("hHGpedRMS","High Gain pedestal RMS",200,0.,20.) ;
270   h31->Sumw2() ;
271   Add2RawsList(h31,kHGpedRMS, !expert, image, !saveCorr) ;
272
273   TH2F * h32 = new TH2F("hpedRMSHGxyMod1","High Gain pedestal RMS Rows x Columns module 1", 64, 0, 64, 56, 0, 56) ;
274   h32->SetOption("colz");
275   Add2RawsList(h32,kHGpedRMSMod1, expert, !image, !saveCorr) ;
276   TH2F * h33 = new TH2F("hpedRMSHGxyMod2","High Gain pedestal RMS Rows x Columns module 2", 64, 0, 64, 56, 0, 56) ;
277   h33->SetOption("colz");
278   Add2RawsList(h33,kHGpedRMSMod2, expert, !image, !saveCorr) ;
279   TH2F * h34 = new TH2F("hpedRMSHGxyMod3","High Gain pedestal RMS Rows x Columns module 3", 64, 0, 64, 56, 0, 56) ;
280   h34->SetOption("colz");
281   Add2RawsList(h34,kHGpedRMSMod3, expert, !image, !saveCorr) ;
282   TH2F * h35 = new TH2F("hpedRMSHGxyMod4","High Gain pedestal RMS Rows x Columns module 4", 64, 0, 64, 56, 0, 56) ;
283   h35->SetOption("colz");
284   Add2RawsList(h35,kHGpedRMSMod4, expert, !image, !saveCorr) ;
285   TH2F * h36 = new TH2F("hpedRMSHGxyMod5","High Gain pedestal RMS Rows x Columns module 5", 64, 0, 64, 56, 0, 56) ;
286   h36->SetOption("colz");
287   Add2RawsList(h36,kHGpedRMSMod5, expert, !image, !saveCorr) ;
288   TH2F * h37 = new TH2F("hpedRMSLGxyMod1","Low Gain pedestal RMS Rows x Columns module 1", 64, 0, 64, 56, 0, 56) ;
289   h37->SetOption("colz");
290   Add2RawsList(h37,kLGpedRMSMod1, expert, !image, !saveCorr) ;
291   TH2F * h38 = new TH2F("hpedRMSLGxyMod2","Low Gain pedestal RMS Rows x Columns module 2", 64, 0, 64, 56, 0, 56) ;
292   h38->SetOption("colz");
293   Add2RawsList(h38,kLGpedRMSMod2, expert, !image, !saveCorr) ;
294   TH2F * h39 = new TH2F("hpedRMSLGxyMod3","Low Gain pedestal RMS Rows x Columns module 3", 64, 0, 64, 56, 0, 56) ;
295   h39->SetOption("colz");
296   Add2RawsList(h39,kLGpedRMSMod3, expert, !image, !saveCorr) ;
297   TH2F * h40 = new TH2F("hpedRMSLGxyMod4","Low Gain pedestal RMS Rows x Columns module 4", 64, 0, 64, 56, 0, 56) ;
298   h40->SetOption("colz");
299   Add2RawsList(h40,kLGpedRMSMod4, expert, !image, !saveCorr) ;
300   TH2F * h41 = new TH2F("hpedRMSLGxyMod5","Low Gain pedestal RMS Rows x Columns module 5", 64, 0, 64, 56, 0, 56) ;
301   h41->SetOption("colz");
302   Add2RawsList(h41,kLGpedRMSMod5, expert, !image, !saveCorr) ;
303
304   /*
305   TH1F * h42 = new TH1F("hLGpedMean","Low Gain pedestal Mean",200,0.,20.) ;
306   h42->Sumw2() ;
307   Add2RawsList(h42,kLGpedMean, !expert, image, !saveCorr) ;
308   TH1F * h43 = new TH1F("hHGpedMean","High Gain pedestal Mean",200,0.,20.) ;
309   h43->Sumw2() ;
310   Add2RawsList(h43,kHGpedMean, !expert, image, !saveCorr) ;
311
312   TH2F * h44 = new TH2F("hpedMeanHGxyMod1","High Gain pedestal Mean Rows x Columns module 1", 64, 0, 64, 56, 0, 56) ;
313   h44->SetOption("colz");
314   Add2RawsList(h44,kHGpedMeanMod1, expert, !image, !saveCorr) ;
315   TH2F * h45 = new TH2F("hpedMeanHGxyMod2","High Gain pedestal Mean Rows x Columns module 2", 64, 0, 64, 56, 0, 56) ;
316   h45->SetOption("colz");
317   Add2RawsList(h45,kHGpedMeanMod2, expert, !image, !saveCorr) ;
318   TH2F * h46 = new TH2F("hpedMeanHGxyMod3","High Gain pedestal Mean Rows x Columns module 3", 64, 0, 64, 56, 0, 56) ;
319   h46->SetOption("colz");
320   Add2RawsList(h46,kHGpedMeanMod3, expert, !image, !saveCorr) ;
321   TH2F * h47 = new TH2F("hpedMeanHGxyMod4","High Gain pedestal Mean Rows x Columns module 4", 64, 0, 64, 56, 0, 56) ;
322   h47->SetOption("colz");
323   Add2RawsList(h47,kHGpedMeanMod4, expert, !image, !saveCorr) ;
324   TH2F * h48 = new TH2F("hpedMeanHGxyMod5","High Gain pedestal Mean Rows x Columns module 5", 64, 0, 64, 56, 0, 56) ;
325   h48->SetOption("colz");
326   Add2RawsList(h48,kHGpedMeanMod5, expert, !image, !saveCorr) ;
327   TH2F * h49 = new TH2F("hpedMeanLGxyMod1","Low Gain pedestal Mean Rows x Columns module 1", 64, 0, 64, 56, 0, 56) ;
328   h49->SetOption("colz");
329   Add2RawsList(h49,kLGpedMeanMod1, expert, !image, !saveCorr) ;
330   TH2F * h50 = new TH2F("hpedMeanLGxyMod2","Low Gain pedestal Mean Rows x Columns module 2", 64, 0, 64, 56, 0, 56) ;
331   h50->SetOption("colz");
332   Add2RawsList(h50,kLGpedMeanMod2, expert, !image, !saveCorr) ;
333   TH2F * h51 = new TH2F("hpedMeanLGxyMod3","Low Gain pedestal Mean Rows x Columns module 3", 64, 0, 64, 56, 0, 56) ;
334   h51->SetOption("colz");
335   Add2RawsList(h51,kLGpedMeanMod3, expert, !image, !saveCorr) ;
336   TH2F * h52 = new TH2F("hpedMeanLGxyMod4","Low Gain pedestal Mean Rows x Columns module 4", 64, 0, 64, 56, 0, 56) ;
337   h52->SetOption("colz");
338   Add2RawsList(h52,kLGpedMeanMod4, expert, !image, !saveCorr) ;
339   TH2F * h53 = new TH2F("hpedMeanLGxyMod5","Low Gain pedestal Mean Rows x Columns module 5", 64, 0, 64, 56, 0, 56) ;
340   h53->SetOption("colz");
341   Add2RawsList(h53,kLGpedMeanMod5, expert, !image, !saveCorr) ;
342   */
343 }
344
345 //____________________________________________________________________________
346 void AliPHOSQADataMakerRec::MakeESDs(AliESDEvent * esd)
347 {
348   // make QA data from ESDs
349
350   Int_t nTot = 0 ; 
351   Double_t eTot = 0 ; 
352   for ( Int_t index = 0; index < esd->GetNumberOfCaloClusters() ; index++ ) {
353     AliESDCaloCluster * clu = esd->GetCaloCluster(index) ;
354     if( clu->IsPHOS() ) {
355       GetESDsData(kESDSpec)->Fill(clu->E()) ;
356       Double_t *pid=clu->GetPid() ;
357       GetESDsData(kESDpid)->Fill(pid[AliPID::kPhoton]) ;
358       eTot+=clu->E() ;
359       nTot++ ;
360     } 
361   }
362   GetESDsData(kESDNtot)->Fill(nTot) ;
363   GetESDsData(kESDEtot)->Fill(eTot) ;
364 }
365
366 //____________________________________________________________________________
367 void AliPHOSQADataMakerRec::MakeRaws(AliRawReader* rawReader)
368 {
369   //Fill prepared histograms with Raw digit properties
370   rawReader->Reset() ;
371   AliPHOSRawDecoder * decoder ;
372   if(strcmp(GetRecoParam()->EMCDecoderVersion(),"v1")==0)
373     decoder=new AliPHOSRawDecoderv1(rawReader);
374   else
375     if(strcmp(GetRecoParam()->EMCDecoderVersion(),"v2")==0)
376       decoder=new AliPHOSRawDecoderv2(rawReader);
377     else
378       decoder=new AliPHOSRawDecoder(rawReader);
379   //decoder->SubtractPedestals(GetRecoParam()->EMCSubtractPedestals());
380   decoder->SubtractPedestals(kTRUE);
381   Double_t lgEtot=0. ;
382   Double_t hgEtot=0. ;
383   Int_t lgNtot=0 ;
384   Int_t hgNtot=0 ;
385
386   while (decoder->NextDigit()) {
387    Int_t module  = decoder->GetModule() ;
388    Int_t row     = decoder->GetRow() ;
389    Int_t col     = decoder->GetColumn() ;
390    Double_t time = decoder->GetTime() ;
391    Double_t energy  = decoder->GetEnergy() ;
392    Bool_t lowGain = decoder->IsLowGain();
393    if (lowGain) {
394      //if(GetRecoParam()->EMCSubtractPedestals())
395        GetRawsData(kLGpedRMS)->Fill(decoder->GetPedestalRMS()) ;
396      if(energy<2.)
397        continue ;
398      switch(module){
399         case 1: GetRawsData(kLGmod1)->Fill(row-0.5,col-0.5) ; break ;
400         case 2: GetRawsData(kLGmod2)->Fill(row-0.5,col-0.5) ; break ;
401         case 3: GetRawsData(kLGmod3)->Fill(row-0.5,col-0.5) ; break ;
402         case 4: GetRawsData(kLGmod4)->Fill(row-0.5,col-0.5) ; break ;
403         case 5: GetRawsData(kLGmod5)->Fill(row-0.5,col-0.5) ; break ;
404      }
405      switch (module){
406         case 1: ((TH2F*)GetRawsData(kLGpedRMSMod1))->Fill(row-0.5,col-0.5,decoder->GetPedestalRMS()) ; break ;
407         case 2: ((TH2F*)GetRawsData(kLGpedRMSMod2))->Fill(row-0.5,col-0.5,decoder->GetPedestalRMS()) ; break ;
408         case 3: ((TH2F*)GetRawsData(kLGpedRMSMod3))->Fill(row-0.5,col-0.5,decoder->GetPedestalRMS()) ; break ;
409         case 4: ((TH2F*)GetRawsData(kLGpedRMSMod4))->Fill(row-0.5,col-0.5,decoder->GetPedestalRMS()) ; break ;
410         case 5: ((TH2F*)GetRawsData(kLGpedRMSMod5))->Fill(row-0.5,col-0.5,decoder->GetPedestalRMS()) ; break ;
411      }
412      //if quality was evaluated, fill histo
413      if(strcmp(GetRecoParam()->EMCDecoderVersion(),"v1")==0){
414        switch (module){
415          case 1: ((TH2F*)GetRawsData(kLGqualMod1))->Fill(row-0.5,col-0.5,decoder->GetSampleQuality()) ; break ;
416          case 2: ((TH2F*)GetRawsData(kLGqualMod2))->Fill(row-0.5,col-0.5,decoder->GetSampleQuality()) ; break ;
417          case 3: ((TH2F*)GetRawsData(kLGqualMod3))->Fill(row-0.5,col-0.5,decoder->GetSampleQuality()) ; break ;
418          case 4: ((TH2F*)GetRawsData(kLGqualMod4))->Fill(row-0.5,col-0.5,decoder->GetSampleQuality()) ; break ;
419          case 5: ((TH2F*)GetRawsData(kLGqualMod5))->Fill(row-0.5,col-0.5,decoder->GetSampleQuality()) ; break ;
420        }
421      }                                  
422      GetRawsData(kNmodLG)->Fill(module) ;
423      GetRawsData(kLGtime)->Fill(time) ; 
424      GetRawsData(kSpecLG)->Fill(energy) ;    
425      lgEtot+=energy ;
426      lgNtot++ ;   
427    } else {        
428      //if this isnon-ZS run - fill pedestal RMS
429      //if(GetRecoParam()->EMCSubtractPedestals())
430        GetRawsData(kHGpedRMS)->Fill(decoder->GetPedestalRMS()) ;
431      if(energy<8.)
432        continue ;
433      switch (module){
434          case 1: GetRawsData(kHGmod1)->Fill(row-0.5,col-0.5) ; break ;
435          case 2: GetRawsData(kHGmod2)->Fill(row-0.5,col-0.5) ; break ;
436          case 3: GetRawsData(kHGmod3)->Fill(row-0.5,col-0.5) ; break ;
437          case 4: GetRawsData(kHGmod4)->Fill(row-0.5,col-0.5) ; break ;
438          case 5: GetRawsData(kHGmod5)->Fill(row-0.5,col-0.5) ; break ;
439      }
440      switch (module){
441          case 1: ((TH2F*)GetRawsData(kHGpedRMSMod1))->Fill(row-0.5,col-0.5,decoder->GetPedestalRMS()) ; break ;
442          case 2: ((TH2F*)GetRawsData(kHGpedRMSMod2))->Fill(row-0.5,col-0.5,decoder->GetPedestalRMS()) ; break ;
443          case 3: ((TH2F*)GetRawsData(kHGpedRMSMod3))->Fill(row-0.5,col-0.5,decoder->GetPedestalRMS()) ; break ;
444          case 4: ((TH2F*)GetRawsData(kHGpedRMSMod4))->Fill(row-0.5,col-0.5,decoder->GetPedestalRMS()) ; break ;
445          case 5: ((TH2F*)GetRawsData(kHGpedRMSMod5))->Fill(row-0.5,col-0.5,decoder->GetPedestalRMS()) ; break ;
446      }               
447      //if quality was evaluated, fill histo
448      if(strcmp(GetRecoParam()->EMCDecoderVersion(),"v1")==0){
449        switch (module){
450          case 1: ((TH2F*)GetRawsData(kHGqualMod1))->Fill(row-0.5,col-0.5,decoder->GetSampleQuality()) ; break ;
451          case 2: ((TH2F*)GetRawsData(kHGqualMod2))->Fill(row-0.5,col-0.5,decoder->GetSampleQuality()) ; break ;
452          case 3: ((TH2F*)GetRawsData(kHGqualMod3))->Fill(row-0.5,col-0.5,decoder->GetSampleQuality()) ; break ;
453          case 4: ((TH2F*)GetRawsData(kHGqualMod4))->Fill(row-0.5,col-0.5,decoder->GetSampleQuality()) ; break ;
454          case 5: ((TH2F*)GetRawsData(kHGqualMod5))->Fill(row-0.5,col-0.5,decoder->GetSampleQuality()) ; break ;
455        }
456      }
457      GetRawsData(kNmodHG)->Fill(module) ; 
458      GetRawsData(kHGtime)->Fill(time) ;  
459      GetRawsData(kSpecHG)->Fill(energy) ;
460      hgEtot+=energy ; 
461      hgNtot++ ;  
462    }                 
463   }                    
464   delete decoder;
465   GetRawsData(kEtotLG)->Fill(lgEtot) ; 
466   TParameter<double> * p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kEtotLG)->GetName()))) ; 
467   if (p) p->SetVal(lgEtot) ; 
468   GetRawsData(kEtotHG)->Fill(hgEtot) ;  
469   p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kEtotHG)->GetName()))) ; 
470   if (p) p->SetVal(hgEtot) ; 
471   GetRawsData(kNtotLG)->Fill(lgNtot) ;
472   p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kNtotLG)->GetName()))) ; 
473   if (p) p->SetVal(lgNtot) ; 
474   GetRawsData(kNtotHG)->Fill(hgNtot) ;
475   p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kNtotHG)->GetName()))) ; 
476   if (p) p->SetVal(hgNtot) ; 
477 }
478
479 //____________________________________________________________________________
480 void AliPHOSQADataMakerRec::MakeDigits(TClonesArray * digits)
481 {
482   // makes data from Digits
483   
484   GetDigitsData(1)->Fill(digits->GetEntriesFast()) ; 
485   TIter next(digits) ; 
486   AliPHOSDigit * digit ; 
487   while ( (digit = dynamic_cast<AliPHOSDigit *>(next())) ) {
488     GetDigitsData(kDigits)->Fill( digit->GetEnergy()) ;
489   }  
490 }
491
492 //____________________________________________________________________________
493 void AliPHOSQADataMakerRec::MakeDigits(TTree * digitTree)
494 {
495         // makes data from Digit Tree
496         TClonesArray * digits = new TClonesArray("AliPHOSDigit", 1000) ; 
497   
498         TBranch * branch = digitTree->GetBranch("PHOS") ;
499         if ( ! branch ) {
500                 AliWarning("PHOS branch in Digit Tree not found") ; 
501         } else {
502                 branch->SetAddress(&digits) ;
503                 branch->GetEntry(0) ; 
504                 MakeDigits(digits) ; 
505         }
506 }
507
508 //____________________________________________________________________________
509 void AliPHOSQADataMakerRec::MakeRecPoints(TTree * clustersTree)
510 {
511   {
512     // makes data from RecPoints
513     TBranch *emcbranch = clustersTree->GetBranch("PHOSEmcRP");
514     if (!emcbranch) { 
515       AliError("can't get the branch with the PHOS EMC clusters !");
516       return;
517     }
518     TObjArray * emcrecpoints = new TObjArray(100) ;
519     emcbranch->SetAddress(&emcrecpoints);
520     emcbranch->GetEntry(0);
521     
522     GetRecPointsData(kRPNtot)->Fill(emcrecpoints->GetEntriesFast()) ; 
523     TIter next(emcrecpoints) ; 
524     AliPHOSEmcRecPoint * rp ; 
525     Double_t eTot = 0. ; 
526     while ( (rp = dynamic_cast<AliPHOSEmcRecPoint *>(next())) ) {
527       GetRecPointsData(kRPSpec)->Fill( rp->GetEnergy()) ;
528       Int_t mod = rp->GetPHOSMod() ;
529       TVector3 pos ;
530       rp->GetLocalPosition(pos) ;
531       switch(mod){
532         case 1: GetRecPointsData(kRPmod1)->Fill(pos.X(),pos.Z()) ; break ;
533         case 2: GetRecPointsData(kRPmod2)->Fill(pos.X(),pos.Z()) ; break ;
534         case 3: GetRecPointsData(kRPmod3)->Fill(pos.X(),pos.Z()) ; break ;
535         case 4: GetRecPointsData(kRPmod4)->Fill(pos.X(),pos.Z()) ; break ;
536         case 5: GetRecPointsData(kRPmod5)->Fill(pos.X(),pos.Z()) ; break ;
537       }
538       eTot+= rp->GetEnergy() ;
539     }
540     GetRecPointsData(kRPEtot)->Fill(eTot) ;
541     emcrecpoints->Delete();
542     delete emcrecpoints;
543   }
544   {
545     TBranch *cpvbranch = clustersTree->GetBranch("PHOSCpvRP");
546     if (!cpvbranch) { 
547       AliError("can't get the branch with the PHOS CPV clusters !");
548       return;
549     }
550     TObjArray *cpvrecpoints = new TObjArray(100) ;
551     cpvbranch->SetAddress(&cpvrecpoints);
552     cpvbranch->GetEntry(0);
553     
554     GetRecPointsData(kRPNcpv)->Fill(cpvrecpoints->GetEntriesFast()) ; 
555     cpvrecpoints->Delete();
556     delete cpvrecpoints;
557   }
558 }
559
560 //____________________________________________________________________________ 
561 void AliPHOSQADataMakerRec::StartOfDetectorCycle()
562 {
563   //Detector specific actions at start of cycle
564   
565 }