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