]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSQADataMakerRec.cxx
printfs removed
[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(kHGqualMod0) && GetRawsData(kHGmod0))
89       GetRawsData(kHGqualMod0)->Divide( GetRawsData(kHGmod0) ) ;
90     if(GetRawsData(kHGqualMod1) && GetRawsData(kHGmod1))
91       GetRawsData(kHGqualMod1)->Divide( GetRawsData(kHGmod1) ) ;
92     if(GetRawsData(kHGqualMod2) && GetRawsData(kHGmod2))
93       GetRawsData(kHGqualMod2)->Divide( GetRawsData(kHGmod2) ) ;
94     if(GetRawsData(kHGqualMod3) && GetRawsData(kHGmod3))
95       GetRawsData(kHGqualMod3)->Divide( GetRawsData(kHGmod3) ) ;
96     if(GetRawsData(kHGqualMod4) && GetRawsData(kHGmod4))
97       GetRawsData(kHGqualMod4)->Divide( GetRawsData(kHGmod4) ) ;
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("hHighPHOSxyMod0","High Gain in PHOS module 0", 64, 0, 64, 56, 0, 56) ;
187   h0->SetXTitle("x, cells"); h0->SetYTitle("z, cells");
188   Add2RawsList(h0,kHGmod0, expert, !image, !saveCorr) ;
189   TH2I * h1 = new TH2I("hHighPHOSxyMod1","High Gain in PHOS module 1", 64, 0, 64, 56, 0, 56) ;
190   h1->SetXTitle("x, cells"); h1->SetYTitle("z, cells");
191   Add2RawsList(h1,kHGmod1, expert, !image, !saveCorr) ;
192   TH2I * h2 = new TH2I("hHighPHOSxyMod2","High Gain in PHOS module 2", 64, 0, 64, 56, 0, 56) ;
193   h2->SetXTitle("x, cells"); h2->SetYTitle("z, cells");
194   Add2RawsList(h2,kHGmod2, expert, !image, !saveCorr) ;
195   TH2I * h3 = new TH2I("hHighPHOSxyMod3","High Gain in PHOS module 3", 64, 0, 64, 56, 0, 56) ;
196   h3->SetXTitle("x, cells"); h3->SetYTitle("z, cells");
197   Add2RawsList(h3,kHGmod3, expert, !image, !saveCorr) ;
198   TH2I * h4 = new TH2I("hHighPHOSxyMod4","High Gain in PHOS module 4", 64, 0, 64, 56, 0, 56) ;
199   h4->SetXTitle("x, cells"); h4->SetYTitle("z, cells");
200   Add2RawsList(h4,kHGmod4, expert, !image, !saveCorr) ;
201
202   TH2I * h5 = new TH2I("hLowPHOSxyMod0","Low Gain in PHOS module 0", 64, 0, 64, 56, 0, 56) ;
203   h5->SetXTitle("x, cells"); h5->SetYTitle("z, cells");
204   Add2RawsList(h5,kLGmod0, expert, !image, !saveCorr) ;
205   TH2I * h6 = new TH2I("hLowPHOSxyMod1","Low Gain in PHOS module 1", 64, 0, 64, 56, 0, 56) ;
206   h6->SetXTitle("x, cells"); h6->SetYTitle("z, cells");
207   Add2RawsList(h6,kLGmod1, expert, !image, !saveCorr) ;
208   TH2I * h7 = new TH2I("hLowPHOSxyMod2","Low Gain in PHOS module 2", 64, 0, 64, 56, 0, 56) ;
209   h7->SetXTitle("x, cells"); h7->SetYTitle("z, cells");
210   Add2RawsList(h7,kLGmod2, expert, !image, !saveCorr) ;
211   TH2I * h8 = new TH2I("hLowPHOSxyMod3","Low Gain in PHOS module 3", 64, 0, 64, 56, 0, 56) ;
212   h8->SetXTitle("x, cells"); h8->SetYTitle("z, cells");
213   Add2RawsList(h8,kLGmod3, expert, !image, !saveCorr) ;
214   TH2I * h9 = new TH2I("hLowPHOSxyMod4","Low Gain in PHOS module 4", 64, 0, 64, 56, 0, 56) ;
215   h9->SetXTitle("x, cells"); h9->SetYTitle("z, cells");
216   Add2RawsList(h9,kLGmod4, expert, !image, !saveCorr) ;
217
218   TH1I * h10 = new TH1I("hLowPhosModules",    "Low Gain Hits in EMCA PHOS modules", 5, 0, 5) ;
219   h10->SetXTitle("Module number");
220   Add2RawsList(h10, kNmodLG, !expert, image, !saveCorr) ;
221   TH1I * h11 = new TH1I("hHighPhosModules",   "High Gain Hits in EMCA PHOS modules",5, 0, 5) ;
222   h11->SetXTitle("Module number");
223   Add2RawsList(h11, kNmodHG, !expert, image, !saveCorr) ;
224
225   TH1F * h12 = new TH1F("hLowPhosRawtime" , "Low Gain Time of raw hits in PHOS" , 500, -50., 200.) ;
226   h12->SetXTitle("Time [samples]");
227   h12->Sumw2() ;
228   Add2RawsList(h12, kLGtime, expert, !image, !saveCorr) ;
229   TH1F * h13 = new TH1F("hHighPhosRawtime", "High Gain Time of raw hits in PHOS", 500, -50., 200.) ;
230   h13->SetXTitle("Time [samples]");
231   h13->Sumw2() ;
232   Add2RawsList(h13, kHGtime, expert, !image, !saveCorr) ;
233
234   TH1F * h14 = new TH1F("hLowPhosRawEnergy" , "Low Gain Energy of raw hits in PHOS" , 512, 0., 1024.) ;
235   h14->SetXTitle("Energy [ADC counts]");
236   h14->Sumw2() ;
237   Add2RawsList(h14, kSpecLG, !expert, image, !saveCorr) ;
238   TH1F * h15 = new TH1F("hHighPhosRawEnergy", "High Gain Energy of raw hits in PHOS", 512, 0., 1024.) ;
239   h15->SetXTitle("Energy [ADC counts]");
240   h15->Sumw2() ;
241   Add2RawsList(h15, kSpecHG, !expert, image, !saveCorr) ;
242
243   TH1F * h16 = new TH1F("hLowNtot" , "Low Gain Total Number of raw hits in PHOS" , 500, 0., 5000.) ;
244   h16->SetXTitle("Number of hits");
245   h16->Sumw2() ;
246   Add2RawsList(h16, kNtotLG, !expert, image, saveCorr) ;
247   TH1F * h17 = new TH1F("hHighNtot", "High Gain Total Number of raw hits in PHOS", 500, 0., 5000.) ;
248   h17->SetXTitle("Number of hits");
249   h17->Sumw2() ;
250   Add2RawsList(h17, kNtotHG, !expert, image, saveCorr) ;
251
252   TH1F * h18 = new TH1F("hLowEtot" , "Low Gain Total Energy of raw hits in PHOS" , 500, 0., 100000.) ;
253   h18->SetXTitle("Energy [ADC counts]");
254   h18->Sumw2() ;
255   Add2RawsList(h18, kEtotLG, !expert, image, saveCorr) ;
256   TH1F * h19 = new TH1F("hHighEtot", "High Gain Total Energy of raw hits in PHOS", 500, 0., 100000.) ;
257   h19->SetXTitle("Energy [ADC counts]");
258   h19->Sumw2() ;
259   Add2RawsList(h19, kEtotHG, !expert, image, saveCorr) ;
260
261   TH2F * h20 = new TH2F("hQualHGxyMod0","High Gain signal quality in module 0", 64, 0, 64, 56, 0, 56) ;
262   h20->SetXTitle("x, cells"); h20->SetYTitle("z, cells");
263   Add2RawsList(h20,kHGqualMod0, expert, !image, !saveCorr) ;
264   TH2F * h21 = new TH2F("hQualHGxyMod1","High Gain signal quality in module 1", 64, 0, 64, 56, 0, 56) ;
265   h21->SetXTitle("x, cells"); h21->SetYTitle("z, cells");
266   Add2RawsList(h21,kHGqualMod1, expert, !image, !saveCorr) ;
267   TH2F * h22 = new TH2F("hQualHGxyMod2","High Gain signal quality in module 2", 64, 0, 64, 56, 0, 56) ;
268   h22->SetXTitle("x, cells"); h22->SetYTitle("z, cells");
269   Add2RawsList(h22,kHGqualMod2, expert, !image, !saveCorr) ;
270   TH2F * h23 = new TH2F("hQualHGxyMod3","High Gain signal quality in module 3", 64, 0, 64, 56, 0, 56) ;
271   h23->SetXTitle("x, cells"); h23->SetYTitle("z, cells");
272   Add2RawsList(h23,kHGqualMod3, expert, !image, !saveCorr) ;
273   TH2F * h24 = new TH2F("hQualHGxyMod4","High Gain signal quality in module 4", 64, 0, 64, 56, 0, 56) ;
274   h24->SetXTitle("x, cells"); h24->SetYTitle("z, cells");
275   Add2RawsList(h24,kHGqualMod4, expert, !image, !saveCorr) ;
276
277   TH2F * h25 = new TH2F("hQualLGxyMod0","Low Gain signal quality in module 0", 64, 0, 64, 56, 0, 56) ;
278   h25->SetXTitle("x, cells"); h25->SetYTitle("z, cells");
279   Add2RawsList(h25,kLGqualMod0, expert, !image, !saveCorr) ;
280   TH2F * h26 = new TH2F("hQualLGxyMod1","Low Gain signal quality in module 1", 64, 0, 64, 56, 0, 56) ;
281   h26->SetXTitle("x, cells"); h26->SetYTitle("z, cells");
282   Add2RawsList(h26,kLGqualMod1, expert, !image, !saveCorr) ;
283   TH2F * h27 = new TH2F("hQualLGxyMod2","Low Gain signal quality in module 2", 64, 0, 64, 56, 0, 56) ;
284   h27->SetXTitle("x, cells"); h27->SetYTitle("z, cells");
285   Add2RawsList(h27,kLGqualMod2, expert, !image, !saveCorr) ;
286   TH2F * h28 = new TH2F("hQualLGxyMod3","Low Gain signal quality in module 3", 64, 0, 64, 56, 0, 56) ;
287   h28->SetXTitle("x, cells"); h28->SetYTitle("z, cells");
288   Add2RawsList(h28,kLGqualMod3, expert, !image, !saveCorr) ;
289   TH2F * h29 = new TH2F("hQualLGxyMod4","Low Gain signal quality in module 4", 64, 0, 64, 56, 0, 56) ;
290   h29->SetXTitle("x, cells"); h29->SetYTitle("z, cells");
291   Add2RawsList(h29,kLGqualMod4, expert, !image, !saveCorr) ;
292
293   TH1F * h30 = new TH1F("hLGpedRMS" ,"Low Gain pedestal RMS" ,200,0.,20.) ;
294   h30->SetXTitle("RMS [ADC counts]");
295   h30->Sumw2() ;
296   Add2RawsList(h30,kLGpedRMS, expert, !image, !saveCorr) ;
297   TH1F * h31 = new TH1F("hHGpedRMS" ,"High Gain pedestal RMS",200,0.,20.) ;
298   h31->SetXTitle("RMS [ADC counts]");
299   h31->Sumw2() ;
300   Add2RawsList(h31,kHGpedRMS, expert, !image, !saveCorr) ;
301
302   TH2F * h32 = new TH2F("hpedRMSHGxyMod0","High Gain pedestal RMS in module 0", 64, 0, 64, 56, 0, 56) ;
303   h32->SetXTitle("x, cells"); h32->SetYTitle("z, cells");
304   Add2RawsList(h32,kHGpedRMSMod0, expert, !image, !saveCorr) ;
305   TH2F * h33 = new TH2F("hpedRMSHGxyMod1","High Gain pedestal RMS in module 1", 64, 0, 64, 56, 0, 56) ;
306   h33->SetXTitle("x, cells"); h33->SetYTitle("z, cells");
307   Add2RawsList(h33,kHGpedRMSMod1, expert, !image, !saveCorr) ;
308   TH2F * h34 = new TH2F("hpedRMSHGxyMod2","High Gain pedestal RMS in module 2", 64, 0, 64, 56, 0, 56) ;
309   h34->SetXTitle("x, cells"); h34->SetYTitle("z, cells");
310   Add2RawsList(h34,kHGpedRMSMod2, expert, !image, !saveCorr) ;
311   TH2F * h35 = new TH2F("hpedRMSHGxyMod3","High Gain pedestal RMS in module 3", 64, 0, 64, 56, 0, 56) ;
312   h35->SetXTitle("x, cells"); h35->SetYTitle("z, cells");
313   Add2RawsList(h35,kHGpedRMSMod3, expert, !image, !saveCorr) ;
314   TH2F * h36 = new TH2F("hpedRMSHGxyMod4","High Gain pedestal RMS in module 4", 64, 0, 64, 56, 0, 56) ;
315   h36->SetXTitle("x, cells"); h36->SetYTitle("z, cells");
316   Add2RawsList(h36,kHGpedRMSMod4, expert, !image, !saveCorr) ;
317
318   TH2F * h37 = new TH2F("hpedRMSLGxyMod0","Low Gain pedestal RMS in module 0", 64, 0, 64, 56, 0, 56) ;
319   h37->SetXTitle("x, cells"); h37->SetYTitle("z, cells");
320   Add2RawsList(h37,kLGpedRMSMod0, expert, !image, !saveCorr) ;
321   TH2F * h38 = new TH2F("hpedRMSLGxyMod1","Low Gain pedestal RMS in module 1", 64, 0, 64, 56, 0, 56) ;
322   h38->SetXTitle("x, cells"); h38->SetYTitle("z, cells");
323   Add2RawsList(h38,kLGpedRMSMod1, expert, !image, !saveCorr) ;
324   TH2F * h39 = new TH2F("hpedRMSLGxyMod2","Low Gain pedestal RMS in module 2", 64, 0, 64, 56, 0, 56) ;
325   h39->SetXTitle("x, cells"); h39->SetYTitle("z, cells");
326   Add2RawsList(h39,kLGpedRMSMod2, expert, !image, !saveCorr) ;
327   TH2F * h40 = new TH2F("hpedRMSLGxyMod3","Low Gain pedestal RMS in module 3", 64, 0, 64, 56, 0, 56) ;
328   h40->SetXTitle("x, cells"); h40->SetYTitle("z, cells");
329   Add2RawsList(h40,kLGpedRMSMod3, expert, !image, !saveCorr) ;
330   TH2F * h41 = new TH2F("hpedRMSLGxyMod4","Low Gain pedestal RMS in module 4", 64, 0, 64, 56, 0, 56) ;
331   h41->SetXTitle("x, cells"); h41->SetYTitle("z, cells");
332   Add2RawsList(h41,kLGpedRMSMod4, expert, !image, !saveCorr) ;
333
334 }
335
336 //____________________________________________________________________________
337 void AliPHOSQADataMakerRec::MakeESDs(AliESDEvent * esd)
338 {
339   // make QA data from ESDs
340  
341   Int_t nTot = 0 ; 
342   Double_t eTot = 0 ; 
343   for ( Int_t index = 0; index < esd->GetNumberOfCaloClusters() ; index++ ) {
344     AliESDCaloCluster * clu = esd->GetCaloCluster(index) ;
345     if( clu->IsPHOS() ) {
346       GetESDsData(kESDSpec)->Fill(clu->E()) ;
347       const Double_t * pid = clu->GetPID() ;
348       GetESDsData(kESDpid)->Fill(pid[AliPID::kPhoton]) ;
349       eTot+=clu->E() ;
350       nTot++ ;
351     } 
352   }
353   GetESDsData(kESDNtot)->Fill(nTot) ;
354   GetESDsData(kESDEtot)->Fill(eTot) ;
355 }
356
357 //____________________________________________________________________________
358 void AliPHOSQADataMakerRec::MakeRaws(AliRawReader* rawReader)
359 {
360   //Fill prepared histograms with Raw digit properties
361
362   rawReader->Reset() ;
363
364   const TObjArray* maps = AliPHOSRecoParam::GetMappings();
365   if(!maps) AliFatal("Cannot retrieve ALTRO mappings!!");
366
367   AliAltroMapping *mapping[20];
368   for(Int_t i = 0; i < 20; i++) {
369     mapping[i] = (AliAltroMapping*)maps->At(i);
370   }
371
372   AliCaloRawStreamV3 fRawStream(rawReader,"PHOS",mapping);
373
374   AliPHOSRawFitterv0 * fitter ;
375   if     (strcmp(GetRecoParam()->EMCFitterVersion(),"v1")==0)
376     fitter=new AliPHOSRawFitterv1();
377   else if(strcmp(GetRecoParam()->EMCFitterVersion(),"v2")==0)
378     fitter=new AliPHOSRawFitterv2();
379   else
380     fitter=new AliPHOSRawFitterv0();
381   Double_t lgEtot=0. ;
382   Double_t hgEtot=0. ;
383   Int_t    lgNtot=0 ;
384   Int_t    hgNtot=0 ;
385
386   while (fRawStream.NextDDL()) {
387     while (fRawStream.NextChannel()) {
388       Int_t module   = fRawStream.GetModule();
389       Int_t cellX    = fRawStream.GetCellX();
390       Int_t cellZ    = fRawStream.GetCellZ();
391       Int_t caloFlag = fRawStream.GetCaloFlag(); // 0=LG, 1=HG, 2=TRU
392
393       if(caloFlag!=0 && caloFlag!=1) continue; //TRU data!
394
395       fitter->SetChannelGeo(module+1,cellX+1,cellZ+1,caloFlag);
396
397       if(fitter->GetAmpOffset()==0 && fitter->GetAmpThreshold()==0) {
398         Short_t altroCFG1 = fRawStream.GetAltroCFG1();
399         Bool_t ZeroSuppressionEnabled = (altroCFG1 >> 15) & 0x1;
400         if(ZeroSuppressionEnabled) {
401           Short_t offset = (altroCFG1 >> 10) & 0xf;
402           Short_t threshold = altroCFG1 & 0x3ff;
403           fitter->SubtractPedestals(kFALSE);
404           fitter->SetAmpOffset(offset);
405           fitter->SetAmpThreshold(threshold);
406         }
407         else
408           fitter->SubtractPedestals(kTRUE);
409       }
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->Eval(sig,sigStart,sigLength);
419       } // End of NextBunch()
420
421       Double_t energy = fitter->GetEnergy() ; 
422       Double_t time   = fitter->GetTime() ;
423
424       if (caloFlag == 0) { // LG
425         GetRawsData(kLGpedRMS)->Fill(fitter->GetPedestalRMS()) ;
426         switch(module){
427         case 0: GetRawsData(kLGmod0)->Fill(cellX,cellZ) ; break ;
428         case 1: GetRawsData(kLGmod1)->Fill(cellX,cellZ) ; break ;
429         case 2: GetRawsData(kLGmod2)->Fill(cellX,cellZ) ; break ;
430         case 3: GetRawsData(kLGmod3)->Fill(cellX,cellZ) ; break ;
431         case 4: GetRawsData(kLGmod4)->Fill(cellX,cellZ) ; break ;
432         }
433         switch (module){
434         case 0: ((TH2F*)GetRawsData(kLGpedRMSMod0))->Fill(cellX,cellZ,fitter->GetPedestalRMS()) ; break ;
435         case 1: ((TH2F*)GetRawsData(kLGpedRMSMod1))->Fill(cellX,cellZ,fitter->GetPedestalRMS()) ; break ;
436         case 2: ((TH2F*)GetRawsData(kLGpedRMSMod2))->Fill(cellX,cellZ,fitter->GetPedestalRMS()) ; break ;
437         case 3: ((TH2F*)GetRawsData(kLGpedRMSMod3))->Fill(cellX,cellZ,fitter->GetPedestalRMS()) ; break ;
438         case 4: ((TH2F*)GetRawsData(kLGpedRMSMod4))->Fill(cellX,cellZ,fitter->GetPedestalRMS()) ; break ;
439         }
440         //if quality was evaluated, fill histo
441         if(strcmp(GetRecoParam()->EMCFitterVersion(),"v1")==0){
442           switch (module){
443           case 0: ((TH2F*)GetRawsData(kLGqualMod0))->Fill(cellX,cellZ,fitter->GetSignalQuality()) ; break ;
444           case 1: ((TH2F*)GetRawsData(kLGqualMod1))->Fill(cellX,cellZ,fitter->GetSignalQuality()) ; break ;
445           case 2: ((TH2F*)GetRawsData(kLGqualMod2))->Fill(cellX,cellZ,fitter->GetSignalQuality()) ; break ;
446           case 3: ((TH2F*)GetRawsData(kLGqualMod3))->Fill(cellX,cellZ,fitter->GetSignalQuality()) ; break ;
447           case 4: ((TH2F*)GetRawsData(kLGqualMod4))->Fill(cellX,cellZ,fitter->GetSignalQuality()) ; break ;
448           }
449         }                                  
450         GetRawsData(kNmodLG)->Fill(module) ;
451         GetRawsData(kLGtime)->Fill(time) ; 
452         GetRawsData(kSpecLG)->Fill(energy) ;    
453         lgEtot+=energy ;
454         lgNtot++ ;   
455       }
456       else if (caloFlag == 1) { // HG        
457         GetRawsData(kHGpedRMS)->Fill(fitter->GetPedestalRMS()) ;
458         switch (module){
459         case 0: GetRawsData(kHGmod0)->Fill(cellX,cellZ) ; break ;
460         case 1: GetRawsData(kHGmod1)->Fill(cellX,cellZ) ; break ;
461         case 2: GetRawsData(kHGmod2)->Fill(cellX,cellZ) ; break ;
462         case 3: GetRawsData(kHGmod3)->Fill(cellX,cellZ) ; break ;
463         case 4: GetRawsData(kHGmod4)->Fill(cellX,cellZ) ; break ;
464         }
465         switch (module){
466         case 0: ((TH2F*)GetRawsData(kHGpedRMSMod0))->Fill(cellX,cellZ,fitter->GetPedestalRMS()) ; break ;
467         case 1: ((TH2F*)GetRawsData(kHGpedRMSMod1))->Fill(cellX,cellZ,fitter->GetPedestalRMS()) ; break ;
468         case 2: ((TH2F*)GetRawsData(kHGpedRMSMod2))->Fill(cellX,cellZ,fitter->GetPedestalRMS()) ; break ;
469         case 3: ((TH2F*)GetRawsData(kHGpedRMSMod3))->Fill(cellX,cellZ,fitter->GetPedestalRMS()) ; break ;
470         case 4: ((TH2F*)GetRawsData(kHGpedRMSMod4))->Fill(cellX,cellZ,fitter->GetPedestalRMS()) ; break ;
471         }               
472         //if quality was evaluated, fill histo
473         if(strcmp(GetRecoParam()->EMCFitterVersion(),"v1")==0){
474           switch (module){
475           case 0: ((TH2F*)GetRawsData(kHGqualMod0))->Fill(cellX,cellZ,fitter->GetSignalQuality()) ; break ;
476           case 1: ((TH2F*)GetRawsData(kHGqualMod1))->Fill(cellX,cellZ,fitter->GetSignalQuality()) ; break ;
477           case 2: ((TH2F*)GetRawsData(kHGqualMod2))->Fill(cellX,cellZ,fitter->GetSignalQuality()) ; break ;
478           case 3: ((TH2F*)GetRawsData(kHGqualMod3))->Fill(cellX,cellZ,fitter->GetSignalQuality()) ; break ;
479           case 4: ((TH2F*)GetRawsData(kHGqualMod4))->Fill(cellX,cellZ,fitter->GetSignalQuality()) ; break ;
480           }       
481         }
482         GetRawsData(kNmodHG)->Fill(module) ; 
483         GetRawsData(kHGtime)->Fill(time) ;  
484         GetRawsData(kSpecHG)->Fill(energy) ;
485         hgEtot+=energy ; 
486         hgNtot++ ;  
487       }
488     }  // End of NextChannel
489   } // End of NextDDL
490   delete fitter;
491
492   GetRawsData(kEtotLG)->Fill(lgEtot) ; 
493   TParameter<double> * p;
494   p = dynamic_cast<TParameter<double>*>(GetParameterList()->
495                                         FindObject(Form("%s_%s_%s", GetName(), 
496                                                         AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), 
497                                                         GetRawsData(kEtotLG)->GetName()))) ; 
498   if (p) p->SetVal(lgEtot) ; 
499   GetRawsData(kEtotHG)->Fill(hgEtot) ;  
500   p = dynamic_cast<TParameter<double>*>(GetParameterList()->
501                                         FindObject(Form("%s_%s_%s", GetName(), 
502                                                         AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), 
503                                                         GetRawsData(kEtotHG)->GetName()))) ; 
504   if (p) p->SetVal(hgEtot) ; 
505   GetRawsData(kNtotLG)->Fill(lgNtot) ;
506   p = dynamic_cast<TParameter<double>*>(GetParameterList()->
507                                         FindObject(Form("%s_%s_%s", GetName(), 
508                                                         AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), 
509                                                         GetRawsData(kNtotLG)->GetName()))) ; 
510   if (p) p->SetVal(lgNtot) ; 
511   GetRawsData(kNtotHG)->Fill(hgNtot) ;
512   p = dynamic_cast<TParameter<double>*>(GetParameterList()->
513                                         FindObject(Form("%s_%s_%s", 
514                                                         GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), 
515                                                         GetRawsData(kNtotHG)->GetName()))) ; 
516   if (p) p->SetVal(hgNtot) ; 
517 }
518
519 //____________________________________________________________________________
520 void AliPHOSQADataMakerRec::MakeDigits()
521 {
522   // makes data from Digits
523
524   TH1 * hist = GetDigitsData(kDigitsMul) ; 
525   if ( ! hist )
526     InitDigits() ;
527   GetDigitsData(kDigitsMul)->Fill(fDigitsArray->GetEntriesFast()) ; 
528   TIter next(fDigitsArray) ; 
529   AliPHOSDigit * digit ; 
530   while ( (digit = dynamic_cast<AliPHOSDigit *>(next())) ) {
531     GetDigitsData(kDigits)->Fill( digit->GetEnergy()) ;
532   }  
533 }
534
535 //____________________________________________________________________________
536 void AliPHOSQADataMakerRec::MakeDigits(TTree * digitTree)
537 {
538         // makes data from Digit Tree
539         if (fDigitsArray) 
540     fDigitsArray->Clear() ; 
541   else 
542     fDigitsArray = new TClonesArray("AliPHOSDigit", 1000) ; 
543   
544         TBranch * branch = digitTree->GetBranch("PHOS") ;
545         if ( ! branch ) {
546                 AliWarning("PHOS branch in Digit Tree not found") ; 
547         } else {
548                 branch->SetAddress(&fDigitsArray) ;
549                 branch->GetEntry(0) ; 
550                 MakeDigits() ; 
551         }
552 }
553
554 //____________________________________________________________________________
555 void AliPHOSQADataMakerRec::MakeRecPoints(TTree * clustersTree)
556 {
557   {
558     // makes data from RecPoints
559     TBranch *emcbranch = clustersTree->GetBranch("PHOSEmcRP");
560     if (!emcbranch) { 
561       AliError("can't get the branch with the PHOS EMC clusters !");
562       return;
563     }
564     
565     TObjArray * emcrecpoints = new TObjArray(100) ;
566     emcbranch->SetAddress(&emcrecpoints);
567     emcbranch->GetEntry(0);
568     
569     GetRecPointsData(kRPNtot)->Fill(emcrecpoints->GetEntriesFast()) ; 
570     TIter next(emcrecpoints) ; 
571     AliPHOSEmcRecPoint * rp ; 
572     Double_t eTot = 0. ; 
573     while ( (rp = static_cast<AliPHOSEmcRecPoint *>(next())) ) {
574       GetRecPointsData(kRPSpec)->Fill( rp->GetEnergy()) ;
575       Int_t mod = rp->GetPHOSMod() ;
576       TVector3 pos ;
577       rp->GetLocalPosition(pos) ;
578       switch(mod){
579         case 1: GetRecPointsData(kRPmod1)->Fill(pos.X(),pos.Z()) ; break ;
580         case 2: GetRecPointsData(kRPmod2)->Fill(pos.X(),pos.Z()) ; break ;
581         case 3: GetRecPointsData(kRPmod3)->Fill(pos.X(),pos.Z()) ; break ;
582         case 4: GetRecPointsData(kRPmod4)->Fill(pos.X(),pos.Z()) ; break ;
583         case 5: GetRecPointsData(kRPmod5)->Fill(pos.X(),pos.Z()) ; break ;
584       }
585       eTot+= rp->GetEnergy() ;
586     }
587     GetRecPointsData(kRPEtot)->Fill(eTot) ;
588     emcrecpoints->Delete();
589     delete emcrecpoints;
590   }
591   {
592     TBranch *cpvbranch = clustersTree->GetBranch("PHOSCpvRP");
593     if (!cpvbranch) { 
594       AliError("can't get the branch with the PHOS CPV clusters !");
595       return;
596     }
597     TObjArray *cpvrecpoints = new TObjArray(100) ;
598     cpvbranch->SetAddress(&cpvrecpoints);
599     cpvbranch->GetEntry(0);
600     
601     GetRecPointsData(kRPNcpv)->Fill(cpvrecpoints->GetEntriesFast()) ; 
602     cpvrecpoints->Delete();
603     delete cpvrecpoints;
604   }
605 }
606
607 //____________________________________________________________________________ 
608 void AliPHOSQADataMakerRec::StartOfDetectorCycle()
609 {
610   //Detector specific actions at start of cycle
611   
612 }