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