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