]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSQADataMakerRec.cxx
add qa check
[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   ResetEventTrigClasses();
85
86   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
87     if (! IsValidEventSpecie(specie, list)) continue;
88     SetEventSpecie(AliRecoParam::ConvertIndex(specie)); 
89     //
90     for (int itc=-1;itc<GetNTrigClasses();itc++) {  // RS: loop over all trigger clones
91       //
92       if(GetRawsData(kHGqualMod0,itc) && GetRawsData(kHGmod0,itc))
93         GetRawsData(kHGqualMod0,itc)->Divide( GetRawsData(kHGmod0,itc) );
94       if(GetRawsData(kHGqualMod1,itc) && GetRawsData(kHGmod1,itc))
95         GetRawsData(kHGqualMod1,itc)->Divide( GetRawsData(kHGmod1,itc) );
96       if(GetRawsData(kHGqualMod2,itc) && GetRawsData(kHGmod2,itc))
97         GetRawsData(kHGqualMod2,itc)->Divide( GetRawsData(kHGmod2,itc) );
98       if(GetRawsData(kHGqualMod3,itc) && GetRawsData(kHGmod3,itc))
99         GetRawsData(kHGqualMod3,itc)->Divide( GetRawsData(kHGmod3,itc) );
100       if(GetRawsData(kHGqualMod4,itc) && GetRawsData(kHGmod4,itc))
101         GetRawsData(kHGqualMod4,itc)->Divide( GetRawsData(kHGmod4,itc) );
102       if(GetRawsData(kNRCUHG,itc)) {
103         Int_t norm = GetEvCountTotalRaws();
104         if (norm < 1) norm = 1;
105         TH1F* tmp = (TH1F*)GetRawsData(kNRCUHGnorm,itc);
106         tmp->Reset();
107         tmp->Add((TH1F*)GetRawsData(kNRCUHG,itc),1./norm);
108       }
109     } // RS: loop over all trigger clones
110   }
111   // do the QA checking
112   AliQAChecker::Instance()->Run(AliQAv1::kPHOS, task, list) ;  
113 }
114
115 //____________________________________________________________________________ 
116 void AliPHOSQADataMakerRec::InitESDs()
117 {
118   //Create histograms to controll ESD
119  
120   const Bool_t expert   = kTRUE ; 
121   const Bool_t image    = kTRUE ; 
122
123   TH1F * h1 = new TH1F("hESDPhosSpectrum",  "ESDs spectrum in PHOS; Energy [MeV];Counts"                ,  200, 0.,   20.) ;
124   h1->Sumw2() ;
125   Add2ESDsList(h1, kESDSpec, !expert, image)  ;
126
127   TH1I * h2 = new TH1I("hESDPhosMul",       "ESDs multiplicity distribution in PHOS; # of clusters;Counts", 100, 0,   100 ) ; 
128   h2->Sumw2() ;
129   Add2ESDsList(h2, kESDNtot, !expert, image) ;
130  
131   TH1F * h3 = new TH1F("hESDPhosEtot",      "ESDs total energy;Energy [MeV];Counts"                     , 2000, 0,  200.) ; 
132   h3->Sumw2() ;
133   Add2ESDsList(h3, kESDEtot, !expert, image) ;  //Expert histo
134  
135   TH1F * h4 = new TH1F("hESDpid",           "ESDs PID distribution in PHOS;Particle Id;Counts"         , 100, 0.,    1.) ;
136   h4->Sumw2() ;
137   Add2ESDsList(h4, kESDpid, !expert, image) ; //Expert histo
138   //
139   ClonePerTrigClass(AliQAv1::kESDS); // this should be the last line    
140 }
141
142 //____________________________________________________________________________ 
143 void AliPHOSQADataMakerRec::InitDigits()
144 {
145   // create Digits histograms in Digits subdir
146   const Bool_t expert   = kTRUE ; 
147   const Bool_t image    = kTRUE ; 
148   TH1I * h0 = new TH1I("hPhosDigits",    "Digits amplitude distribution in PHOS;Amplitude [ADC counts];Counts",    500, 0, 1000) ; 
149   h0->Sumw2() ;
150   Add2DigitsList(h0, kDigits, !expert, image) ;
151   TH1I * h1 = new TH1I("hPhosDigitsMul", "Digits multiplicity distribution in PHOS;# of Digits;Entries", 2000, 0, 10000) ; 
152   h1->Sumw2() ;
153   Add2DigitsList(h1, kDigitsMul, !expert, image) ;
154   //
155   ClonePerTrigClass(AliQAv1::kDIGITS); // this should be the last line
156 }
157
158 //____________________________________________________________________________ 
159 void AliPHOSQADataMakerRec::InitRecPoints()
160 {
161   // create Reconstructed Points histograms in RecPoints subdir
162   const Bool_t expert   = kTRUE ; 
163   const Bool_t image    = kTRUE ; 
164   
165   TH2I * h0 = new TH2I("hRpPHOSxyMod1","RecPoints Rows x Columns for PHOS module 1;Row #;Column #", 64, -72., 72., 56, -63., 63.) ;                             
166   Add2RecPointsList(h0,kRPmod1, expert, !image) ;
167   TH2I * h1 = new TH2I("hRpPHOSxyMod2","RecPoints Rows x Columns for PHOS module 2;Row #;Column #", 64, -72., 72., 56, -63., 63.) ;                             
168   Add2RecPointsList(h1,kRPmod2, expert, !image) ;
169   TH2I * h2 = new TH2I("hRpPHOSxyMod3","RecPoints Rows x Columns for PHOS module 3;Row #;Column #", 64, -72., 72., 56, -63., 63.) ;                             
170   Add2RecPointsList(h2,kRPmod3, expert, !image) ;
171   TH2I * h3 = new TH2I("hRpPHOSxyMod4","RecPoints Rows x Columns for PHOS module 4;Row #;Column #", 64, -72., 72., 56, -63., 63.) ;                             
172   Add2RecPointsList(h3,kRPmod4, expert, !image) ;
173   TH2I * h4 = new TH2I("hRpPHOSxyMod5","RecPoints Rows x Columns for PHOS module 5;Row #;Column #", 64, -72., 72., 56, -63., 63.) ;                             
174   Add2RecPointsList(h4,kRPmod5, expert, !image) ;
175  
176   TH1F * h5 = new TH1F("hEmcPhosRecPointsSpectrum",  "EMC RecPoints spectrum in PHOS;Energy [MeV];Counts",   2000, 0., 20.) ; 
177   h5->Sumw2() ;
178   Add2RecPointsList(h5, kRPSpec, !expert, image)  ;
179
180   TH1I * h6 = new TH1I("hEmcPhosRecPointsMul", "EMC RecPoints multiplicity distribution in PHOS;# of EMC Clusters;Entries", 100, 0,  100) ; 
181   h6->Sumw2() ;
182   Add2RecPointsList(h6, kRPNtot, !expert, image) ;
183
184   TH1I * h7 = new TH1I("hEmcPhosRecPointsEtot", "EMC RecPoints Etot;Energy [MeV];Counts", 200, 0,  200.) ; 
185   h7->Sumw2() ;
186   Add2RecPointsList(h7, kRPEtot, !expert, image) ;
187
188   TH1I * h8 = new TH1I("hCpvPhosRecPointsMul", "CPV RecPoints multiplicity distribution in PHOS;# of CPV clusters;Counts", 100, 0,  100) ; 
189   h8->Sumw2() ;
190   Add2RecPointsList(h8, kRPNcpv, !expert, image) ;
191   //
192   ClonePerTrigClass(AliQAv1::kRECPOINTS); // this should be the last line
193 }
194
195 //____________________________________________________________________________ 
196 void AliPHOSQADataMakerRec::InitRaws()
197 {
198   // create Raws histograms in Raws subdir
199   const Bool_t expert   = kTRUE ; 
200   const Bool_t saveCorr = kTRUE ; 
201   const Bool_t image    = kTRUE ; 
202
203   TH2I * h0 = new TH2I("hHighPHOSxyMod0","High Gain in PHOS module 0", 64, 0, 64, 56, 0, 56) ;
204   h0->SetXTitle("x, cells"); h0->SetYTitle("z, cells");
205   Add2RawsList(h0,kHGmod0, expert, !image, !saveCorr) ;
206   TH2I * h1 = new TH2I("hHighPHOSxyMod1","High Gain in PHOS module 1", 64, 0, 64, 56, 0, 56) ;
207   h1->SetXTitle("x, cells"); h1->SetYTitle("z, cells");
208   Add2RawsList(h1,kHGmod1, expert, !image, !saveCorr) ;
209   TH2I * h2 = new TH2I("hHighPHOSxyMod2","High Gain in PHOS module 2", 64, 0, 64, 56, 0, 56) ;
210   h2->SetXTitle("x, cells"); h2->SetYTitle("z, cells");
211   Add2RawsList(h2,kHGmod2, expert, !image, !saveCorr) ;
212   TH2I * h3 = new TH2I("hHighPHOSxyMod3","High Gain in PHOS module 3", 64, 0, 64, 56, 0, 56) ;
213   h3->SetXTitle("x, cells"); h3->SetYTitle("z, cells");
214   Add2RawsList(h3,kHGmod3, expert, !image, !saveCorr) ;
215   TH2I * h4 = new TH2I("hHighPHOSxyMod4","High Gain in PHOS module 4", 64, 0, 64, 56, 0, 56) ;
216   h4->SetXTitle("x, cells"); h4->SetYTitle("z, cells");
217   Add2RawsList(h4,kHGmod4, expert, !image, !saveCorr) ;
218
219   TH2I * h5 = new TH2I("hLowPHOSxyMod0","Low Gain in PHOS module 0", 64, 0, 64, 56, 0, 56) ;
220   h5->SetXTitle("x, cells"); h5->SetYTitle("z, cells");
221   Add2RawsList(h5,kLGmod0, expert, !image, !saveCorr) ;
222   TH2I * h6 = new TH2I("hLowPHOSxyMod1","Low Gain in PHOS module 1", 64, 0, 64, 56, 0, 56) ;
223   h6->SetXTitle("x, cells"); h6->SetYTitle("z, cells");
224   Add2RawsList(h6,kLGmod1, expert, !image, !saveCorr) ;
225   TH2I * h7 = new TH2I("hLowPHOSxyMod2","Low Gain in PHOS module 2", 64, 0, 64, 56, 0, 56) ;
226   h7->SetXTitle("x, cells"); h7->SetYTitle("z, cells");
227   Add2RawsList(h7,kLGmod2, expert, !image, !saveCorr) ;
228   TH2I * h8 = new TH2I("hLowPHOSxyMod3","Low Gain in PHOS module 3", 64, 0, 64, 56, 0, 56) ;
229   h8->SetXTitle("x, cells"); h8->SetYTitle("z, cells");
230   Add2RawsList(h8,kLGmod3, expert, !image, !saveCorr) ;
231   TH2I * h9 = new TH2I("hLowPHOSxyMod4","Low Gain in PHOS module 4", 64, 0, 64, 56, 0, 56) ;
232   h9->SetXTitle("x, cells"); h9->SetYTitle("z, cells");
233   Add2RawsList(h9,kLGmod4, expert, !image, !saveCorr) ;
234
235   TH1I * h10 = new TH1I("hLowPhosModules",    "Low Gain Hits in EMCA PHOS modules", 5, 0, 5) ;
236   h10->SetXTitle("Module number");
237   Add2RawsList(h10, kNmodLG, !expert, image, !saveCorr) ;
238   TH1I * h11 = new TH1I("hHighPhosModules",   "High Gain Hits in EMCA PHOS modules",5, 0, 5) ;
239   h11->SetXTitle("Module number");
240   Add2RawsList(h11, kNmodHG, !expert, image, !saveCorr) ;
241
242   TH1I * h11_RCU = new TH1I("hHighPhosRCU",   "PHOS RCU occupancy with HG channels",12, 0, 12) ;
243   h11_RCU->SetXTitle("RCU number");
244   Add2RawsList(h11_RCU, kNRCUHG, expert, image, !saveCorr) ;
245   TH1F * h11_RCUnorm = new TH1F("hHighPhosRCUnorm",   "PHOS RCU normalized occupancy with HG channels",12, 0, 12) ;
246   h11_RCUnorm->SetXTitle("RCU number");
247   h11_RCUnorm->SetYTitle("N hits per event");
248   h11_RCUnorm->Sumw2();
249   h11_RCUnorm->SetMarkerStyle(20);
250   Add2RawsList(h11_RCUnorm, kNRCUHGnorm, !expert, image, !saveCorr) ;
251
252   TH1F * h12 = new TH1F("hLowPhosRawtime" , "Low Gain Time of raw hits in PHOS" , 500, -50., 200.) ;
253   h12->SetXTitle("Time [samples]");
254   h12->Sumw2() ;
255   Add2RawsList(h12, kLGtime, expert, !image, !saveCorr) ;
256   TH1F * h13 = new TH1F("hHighPhosRawtime", "High Gain Time of raw hits in PHOS", 500, -50., 200.) ;
257   h13->SetXTitle("Time [samples]");
258   h13->Sumw2() ;
259   Add2RawsList(h13, kHGtime, expert, !image, !saveCorr) ;
260
261   TH1F * h14 = new TH1F("hLowPhosRawEnergy" , "Low Gain Energy of raw hits in PHOS" , 512, 0., 1024.) ;
262   h14->SetXTitle("Energy [ADC counts]");
263   h14->Sumw2() ;
264   Add2RawsList(h14, kSpecLG, !expert, image, !saveCorr) ;
265   TH1F * h15 = new TH1F("hHighPhosRawEnergy", "High Gain Energy of raw hits in PHOS", 512, 0., 1024.) ;
266   h15->SetXTitle("Energy [ADC counts]");
267   h15->Sumw2() ;
268   Add2RawsList(h15, kSpecHG, !expert, image, !saveCorr) ;
269
270   TH1F * h16 = new TH1F("hLowNtot" , "Low Gain Total Number of raw hits in PHOS" , 500, 0., 5000.) ;
271   h16->SetXTitle("Number of hits");
272   h16->Sumw2() ;
273   Add2RawsList(h16, kNtotLG, !expert, image, saveCorr) ;
274   TH1F * h17 = new TH1F("hHighNtot", "High Gain Total Number of raw hits in PHOS", 500, 0., 5000.) ;
275   h17->SetXTitle("Number of hits");
276   h17->Sumw2() ;
277   Add2RawsList(h17, kNtotHG, !expert, image, saveCorr) ;
278
279   TH1F * h17_1_0 = new TH1F("hHighNtot_1_0", "High Gain Total Number of raw hits in module 1 RCU0", 500, 0., 1000.) ;
280   h17_1_0->SetXTitle("Number of hits");
281   h17_1_0->Sumw2() ;
282   Add2RawsList(h17_1_0, kNtotHG_1_0, expert, image, saveCorr) ;
283
284   TH1F * h17_1_1 = new TH1F("hHighNtot_1_1", "High Gain Total Number of raw hits in module 1 RCU1", 500, 0., 1000.) ;
285   h17_1_1->SetXTitle("Number of hits");
286   h17_1_1->Sumw2() ;
287   Add2RawsList(h17_1_1, kNtotHG_1_1, expert, image, saveCorr) ;
288
289   TH1F * h17_1_2 = new TH1F("hHighNtot_1_2", "High Gain Total Number of raw hits in module 1 RCU2", 500, 0., 1000.) ;
290   h17_1_2->SetXTitle("Number of hits");
291   h17_1_2->Sumw2() ;
292   Add2RawsList(h17_1_2, kNtotHG_1_2, expert, image, saveCorr) ;
293
294   TH1F * h17_1_3 = new TH1F("hHighNtot_1_3", "High Gain Total Number of raw hits in module 1 RCU3", 500, 0., 1000.) ;
295   h17_1_3->SetXTitle("Number of hits");
296   h17_1_3->Sumw2() ;
297   Add2RawsList(h17_1_3, kNtotHG_1_3, expert, image, saveCorr) ;
298
299   TH1F * h17_2_0 = new TH1F("hHighNtot_2_0", "High Gain Total Number of raw hits in module 2 RCU0", 500, 0., 1000.) ;
300   h17_2_0->SetXTitle("Number of hits");
301   h17_2_0->Sumw2() ;
302   Add2RawsList(h17_2_0, kNtotHG_2_0, expert, image, saveCorr) ;
303
304   TH1F * h17_2_1 = new TH1F("hHighNtot_2_1", "High Gain Total Number of raw hits in module 2 RCU1", 500, 0., 1000.) ;
305   h17_2_1->SetXTitle("Number of hits");
306   h17_2_1->Sumw2() ;
307   Add2RawsList(h17_2_1, kNtotHG_2_1, expert, image, saveCorr) ;
308
309   TH1F * h17_2_2 = new TH1F("hHighNtot_2_2", "High Gain Total Number of raw hits in module 2 RCU2", 500, 0., 1000.) ;
310   h17_2_2->SetXTitle("Number of hits");
311   h17_2_2->Sumw2() ;
312   Add2RawsList(h17_2_2, kNtotHG_2_2, expert, image, saveCorr) ;
313
314   TH1F * h17_2_3 = new TH1F("hHighNtot_2_3", "High Gain Total Number of raw hits in module 2 RCU3", 500, 0., 1000.) ;
315   h17_2_3->SetXTitle("Number of hits");
316   h17_2_3->Sumw2() ;
317   Add2RawsList(h17_2_3, kNtotHG_2_3, expert, image, saveCorr) ;
318
319   TH1F * h17_3_0 = new TH1F("hHighNtot_3_0", "High Gain Total Number of raw hits in module 3 RCU0", 500, 0., 1000.) ;
320   h17_3_0->SetXTitle("Number of hits");
321   h17_3_0->Sumw2() ;
322   Add2RawsList(h17_3_0, kNtotHG_3_0, expert, image, saveCorr) ;
323
324   TH1F * h17_3_1 = new TH1F("hHighNtot_3_1", "High Gain Total Number of raw hits in module 3 RCU1", 500, 0., 1000.) ;
325   h17_3_1->SetXTitle("Number of hits");
326   h17_3_1->Sumw2() ;
327   Add2RawsList(h17_3_1, kNtotHG_3_1, expert, image, saveCorr) ;
328
329   TH1F * h17_3_2 = new TH1F("hHighNtot_3_2", "High Gain Total Number of raw hits in module 3 RCU2", 500, 0., 1000.) ;
330   h17_3_2->SetXTitle("Number of hits");
331   h17_3_2->Sumw2() ;
332   Add2RawsList(h17_3_2, kNtotHG_3_2, expert, image, saveCorr) ;
333
334   TH1F * h17_3_3 = new TH1F("hHighNtot_3_3", "High Gain Total Number of raw hits in module 3 RCU3", 500, 0., 1000.) ;
335   h17_3_3->SetXTitle("Number of hits");
336   h17_3_3->Sumw2() ;
337   Add2RawsList(h17_3_3, kNtotHG_3_3, expert, image, saveCorr) ;
338
339   TH1F * h18 = new TH1F("hLowEtot" , "Low Gain Total Energy of raw hits in PHOS" , 500, 0., 100000.) ;
340   h18->SetXTitle("Energy [ADC counts]");
341   h18->Sumw2() ;
342   Add2RawsList(h18, kEtotLG, !expert, image, saveCorr) ;
343   TH1F * h19 = new TH1F("hHighEtot", "High Gain Total Energy of raw hits in PHOS", 500, 0., 100000.) ;
344   h19->SetXTitle("Energy [ADC counts]");
345   h19->Sumw2() ;
346   Add2RawsList(h19, kEtotHG, !expert, image, saveCorr) ;
347
348   TH2F * h20 = new TH2F("hQualHGxyMod0","High Gain signal quality in module 0", 64, 0, 64, 56, 0, 56) ;
349   h20->SetXTitle("x, cells"); h20->SetYTitle("z, cells");
350   Add2RawsList(h20,kHGqualMod0, expert, !image, !saveCorr) ;
351   TH2F * h21 = new TH2F("hQualHGxyMod1","High Gain signal quality in module 1", 64, 0, 64, 56, 0, 56) ;
352   h21->SetXTitle("x, cells"); h21->SetYTitle("z, cells");
353   Add2RawsList(h21,kHGqualMod1, expert, !image, !saveCorr) ;
354   TH2F * h22 = new TH2F("hQualHGxyMod2","High Gain signal quality in module 2", 64, 0, 64, 56, 0, 56) ;
355   h22->SetXTitle("x, cells"); h22->SetYTitle("z, cells");
356   Add2RawsList(h22,kHGqualMod2, expert, !image, !saveCorr) ;
357   TH2F * h23 = new TH2F("hQualHGxyMod3","High Gain signal quality in module 3", 64, 0, 64, 56, 0, 56) ;
358   h23->SetXTitle("x, cells"); h23->SetYTitle("z, cells");
359   Add2RawsList(h23,kHGqualMod3, expert, !image, !saveCorr) ;
360   TH2F * h24 = new TH2F("hQualHGxyMod4","High Gain signal quality in module 4", 64, 0, 64, 56, 0, 56) ;
361   h24->SetXTitle("x, cells"); h24->SetYTitle("z, cells");
362   Add2RawsList(h24,kHGqualMod4, expert, !image, !saveCorr) ;
363
364   TH2F * h25 = new TH2F("hQualLGxyMod0","Low Gain signal quality in module 0", 64, 0, 64, 56, 0, 56) ;
365   h25->SetXTitle("x, cells"); h25->SetYTitle("z, cells");
366   Add2RawsList(h25,kLGqualMod0, expert, !image, !saveCorr) ;
367   TH2F * h26 = new TH2F("hQualLGxyMod1","Low Gain signal quality in module 1", 64, 0, 64, 56, 0, 56) ;
368   h26->SetXTitle("x, cells"); h26->SetYTitle("z, cells");
369   Add2RawsList(h26,kLGqualMod1, expert, !image, !saveCorr) ;
370   TH2F * h27 = new TH2F("hQualLGxyMod2","Low Gain signal quality in module 2", 64, 0, 64, 56, 0, 56) ;
371   h27->SetXTitle("x, cells"); h27->SetYTitle("z, cells");
372   Add2RawsList(h27,kLGqualMod2, expert, !image, !saveCorr) ;
373   TH2F * h28 = new TH2F("hQualLGxyMod3","Low Gain signal quality in module 3", 64, 0, 64, 56, 0, 56) ;
374   h28->SetXTitle("x, cells"); h28->SetYTitle("z, cells");
375   Add2RawsList(h28,kLGqualMod3, expert, !image, !saveCorr) ;
376   TH2F * h29 = new TH2F("hQualLGxyMod4","Low Gain signal quality in module 4", 64, 0, 64, 56, 0, 56) ;
377   h29->SetXTitle("x, cells"); h29->SetYTitle("z, cells");
378   Add2RawsList(h29,kLGqualMod4, expert, !image, !saveCorr) ;
379
380   TH1F * h30 = new TH1F("hLGpedRMS" ,"Low Gain pedestal RMS" ,200,0.,20.) ;
381   h30->SetXTitle("RMS [ADC counts]");
382   h30->Sumw2() ;
383   Add2RawsList(h30,kLGpedRMS, expert, !image, !saveCorr) ;
384   TH1F * h31 = new TH1F("hHGpedRMS" ,"High Gain pedestal RMS",200,0.,20.) ;
385   h31->SetXTitle("RMS [ADC counts]");
386   h31->Sumw2() ;
387   Add2RawsList(h31,kHGpedRMS, expert, !image, !saveCorr) ;
388
389   TH2F * h32 = new TH2F("hpedRMSHGxyMod0","High Gain pedestal RMS in module 0", 64, 0, 64, 56, 0, 56) ;
390   h32->SetXTitle("x, cells"); h32->SetYTitle("z, cells");
391   Add2RawsList(h32,kHGpedRMSMod0, expert, !image, !saveCorr) ;
392   TH2F * h33 = new TH2F("hpedRMSHGxyMod1","High Gain pedestal RMS in module 1", 64, 0, 64, 56, 0, 56) ;
393   h33->SetXTitle("x, cells"); h33->SetYTitle("z, cells");
394   Add2RawsList(h33,kHGpedRMSMod1, expert, !image, !saveCorr) ;
395   TH2F * h34 = new TH2F("hpedRMSHGxyMod2","High Gain pedestal RMS in module 2", 64, 0, 64, 56, 0, 56) ;
396   h34->SetXTitle("x, cells"); h34->SetYTitle("z, cells");
397   Add2RawsList(h34,kHGpedRMSMod2, expert, !image, !saveCorr) ;
398   TH2F * h35 = new TH2F("hpedRMSHGxyMod3","High Gain pedestal RMS in module 3", 64, 0, 64, 56, 0, 56) ;
399   h35->SetXTitle("x, cells"); h35->SetYTitle("z, cells");
400   Add2RawsList(h35,kHGpedRMSMod3, expert, !image, !saveCorr) ;
401   TH2F * h36 = new TH2F("hpedRMSHGxyMod4","High Gain pedestal RMS in module 4", 64, 0, 64, 56, 0, 56) ;
402   h36->SetXTitle("x, cells"); h36->SetYTitle("z, cells");
403   Add2RawsList(h36,kHGpedRMSMod4, expert, !image, !saveCorr) ;
404
405   TH2F * h37 = new TH2F("hpedRMSLGxyMod0","Low Gain pedestal RMS in module 0", 64, 0, 64, 56, 0, 56) ;
406   h37->SetXTitle("x, cells"); h37->SetYTitle("z, cells");
407   Add2RawsList(h37,kLGpedRMSMod0, expert, !image, !saveCorr) ;
408   TH2F * h38 = new TH2F("hpedRMSLGxyMod1","Low Gain pedestal RMS in module 1", 64, 0, 64, 56, 0, 56) ;
409   h38->SetXTitle("x, cells"); h38->SetYTitle("z, cells");
410   Add2RawsList(h38,kLGpedRMSMod1, expert, !image, !saveCorr) ;
411   TH2F * h39 = new TH2F("hpedRMSLGxyMod2","Low Gain pedestal RMS in module 2", 64, 0, 64, 56, 0, 56) ;
412   h39->SetXTitle("x, cells"); h39->SetYTitle("z, cells");
413   Add2RawsList(h39,kLGpedRMSMod2, expert, !image, !saveCorr) ;
414   TH2F * h40 = new TH2F("hpedRMSLGxyMod3","Low Gain pedestal RMS in module 3", 64, 0, 64, 56, 0, 56) ;
415   h40->SetXTitle("x, cells"); h40->SetYTitle("z, cells");
416   Add2RawsList(h40,kLGpedRMSMod3, expert, !image, !saveCorr) ;
417   TH2F * h41 = new TH2F("hpedRMSLGxyMod4","Low Gain pedestal RMS in module 4", 64, 0, 64, 56, 0, 56) ;
418   h41->SetXTitle("x, cells"); h41->SetYTitle("z, cells");
419   Add2RawsList(h41,kLGpedRMSMod4, expert, !image, !saveCorr) ;
420   //
421   ClonePerTrigClass(AliQAv1::kRAWS); // this should be the last line
422 }
423
424 //____________________________________________________________________________
425 void AliPHOSQADataMakerRec::MakeESDs(AliESDEvent * esd)
426 {
427   // make QA data from ESDs 
428   Int_t nTot = 0 ; 
429   Double_t eTot = 0 ; 
430   for ( Int_t index = 0; index < esd->GetNumberOfCaloClusters() ; index++ ) {
431     AliESDCaloCluster * clu = esd->GetCaloCluster(index) ;
432     if( clu->IsPHOS() ) {
433       FillESDsData(kESDSpec,clu->E()) ;
434       const Double_t * pid = clu->GetPID() ;
435       FillESDsData(kESDpid,pid[AliPID::kPhoton]) ;
436       eTot+=clu->E() ;
437       nTot++ ;
438     } 
439   }
440   FillESDsData(kESDNtot,nTot) ;
441   FillESDsData(kESDEtot,eTot) ;
442   //
443   IncEvCountCycleESDs();
444   IncEvCountTotalESDs();
445   //
446 }
447
448 //____________________________________________________________________________
449 void AliPHOSQADataMakerRec::MakeRaws(AliRawReader* rawReader)
450 {
451   //Fill prepared histograms with Raw digit properties
452
453   rawReader->Reset() ;
454
455   const TObjArray* maps = AliPHOSRecoParam::GetMappings();
456   if(!maps) AliFatal("Cannot retrieve ALTRO mappings!!");
457
458   AliAltroMapping *mapping[20];
459   for(Int_t i = 0; i < 20; i++) {
460     mapping[i] = (AliAltroMapping*)maps->At(i);
461   }
462
463   AliCaloRawStreamV3 fRawStream(rawReader,"PHOS",mapping);
464
465   AliPHOSRawFitterv0 * fitter ;
466   if     (strcmp(GetRecoParam()->EMCFitterVersion(),"v1")==0)
467     fitter=new AliPHOSRawFitterv1();
468   else if(strcmp(GetRecoParam()->EMCFitterVersion(),"v2")==0)
469     fitter=new AliPHOSRawFitterv2();
470   else
471     fitter=new AliPHOSRawFitterv0();
472   Double_t lgEtot=0. ;
473   Double_t hgEtot=0. ;
474   Int_t    lgNtot=0 ;
475   Int_t    hgNtot=0 ;
476   Int_t    hgNtotM1RCU0=0, hgNtotM1RCU1=0, hgNtotM1RCU2=0, hgNtotM1RCU3=0;
477   Int_t    hgNtotM2RCU0=0, hgNtotM2RCU1=0, hgNtotM2RCU2=0, hgNtotM2RCU3=0;
478   Int_t    hgNtotM3RCU0=0, hgNtotM3RCU1=0, hgNtotM3RCU2=0, hgNtotM3RCU3=0;
479
480
481   while (fRawStream.NextDDL()) {
482     Int_t RCUnum = fRawStream.GetDDLNumber() - 8;
483
484     while (fRawStream.NextChannel()) {
485       Int_t module   = fRawStream.GetModule();
486       Int_t cellX    = fRawStream.GetCellX();
487       Int_t cellZ    = fRawStream.GetCellZ();
488       Int_t caloFlag = fRawStream.GetCaloFlag(); // 0=LG, 1=HG, 2=TRU
489
490       if(caloFlag!=0 && caloFlag!=1) continue; //TRU data!
491
492       fitter->SetChannelGeo(module+1,cellX+1,cellZ+1,caloFlag);
493
494       if(fitter->GetAmpOffset()==0 && fitter->GetAmpThreshold()==0) {
495         Short_t altroCFG1 = fRawStream.GetAltroCFG1();
496         Bool_t ZeroSuppressionEnabled = (altroCFG1 >> 15) & 0x1;
497         if(ZeroSuppressionEnabled) {
498           Short_t offset = (altroCFG1 >> 10) & 0xf;
499           Short_t threshold = altroCFG1 & 0x3ff;
500           fitter->SubtractPedestals(kFALSE);
501           fitter->SetAmpOffset(offset);
502           fitter->SetAmpThreshold(threshold);
503         }
504         else
505           fitter->SubtractPedestals(kTRUE);
506       }
507
508       Int_t nBunches = 0;
509       while (fRawStream.NextBunch()) {
510         nBunches++;
511         if (nBunches > 1) continue;
512         const UShort_t *sig = fRawStream.GetSignals();
513         Int_t sigStart      = fRawStream.GetStartTimeBin();
514         Int_t sigLength     = fRawStream.GetBunchLength();
515         fitter->Eval(sig,sigStart,sigLength);
516       } // End of NextBunch()
517
518       Double_t energy = fitter->GetEnergy() ; 
519       Double_t time   = fitter->GetTime() ;
520
521       if (caloFlag == 0) { // LG
522         FillRawsData(kLGpedRMS,fitter->GetPedestalRMS()) ;
523         switch(module){
524         case 0: FillRawsData(kLGmod0,cellX,cellZ) ; break ;
525         case 1: FillRawsData(kLGmod1,cellX,cellZ) ; break ;
526         case 2: FillRawsData(kLGmod2,cellX,cellZ) ; break ;
527         case 3: FillRawsData(kLGmod3,cellX,cellZ) ; break ;
528         case 4: FillRawsData(kLGmod4,cellX,cellZ) ; break ;
529         }
530         switch (module){
531         case 0: FillRawsData(kLGpedRMSMod0,cellX,cellZ,fitter->GetPedestalRMS()) ; break ;
532         case 1: FillRawsData(kLGpedRMSMod1,cellX,cellZ,fitter->GetPedestalRMS()) ; break ;
533         case 2: FillRawsData(kLGpedRMSMod2,cellX,cellZ,fitter->GetPedestalRMS()) ; break ;
534         case 3: FillRawsData(kLGpedRMSMod3,cellX,cellZ,fitter->GetPedestalRMS()) ; break ;
535         case 4: FillRawsData(kLGpedRMSMod4,cellX,cellZ,fitter->GetPedestalRMS()) ; break ;
536         }
537         //if quality was evaluated, fill histo
538         if(strcmp(GetRecoParam()->EMCFitterVersion(),"v1")==0){
539           switch (module){
540           case 0: FillRawsData(kLGqualMod0,cellX,cellZ,fitter->GetSignalQuality()) ; break ;
541           case 1: FillRawsData(kLGqualMod1,cellX,cellZ,fitter->GetSignalQuality()) ; break ;
542           case 2: FillRawsData(kLGqualMod2,cellX,cellZ,fitter->GetSignalQuality()) ; break ;
543           case 3: FillRawsData(kLGqualMod3,cellX,cellZ,fitter->GetSignalQuality()) ; break ;
544           case 4: FillRawsData(kLGqualMod4,cellX,cellZ,fitter->GetSignalQuality()) ; break ;
545           }
546         }                                  
547         FillRawsData(kNmodLG,module) ;
548         FillRawsData(kLGtime,time) ; 
549         FillRawsData(kSpecLG,energy) ;    
550         lgEtot+=energy ;
551         lgNtot++ ;   
552       }
553       else if (caloFlag == 1) { // HG        
554         FillRawsData(kHGpedRMS,fitter->GetPedestalRMS()) ;
555         switch (module){
556         case 0: FillRawsData(kHGmod0,cellX,cellZ) ; break ;
557         case 1: FillRawsData(kHGmod1,cellX,cellZ) ; break ;
558         case 2: FillRawsData(kHGmod2,cellX,cellZ) ; break ;
559         case 3: FillRawsData(kHGmod3,cellX,cellZ) ; break ;
560         case 4: FillRawsData(kHGmod4,cellX,cellZ) ; break ;
561         }
562         switch (module){
563         case 0: FillRawsData(kHGpedRMSMod0,cellX,cellZ,fitter->GetPedestalRMS()) ; break ;
564         case 1: FillRawsData(kHGpedRMSMod1,cellX,cellZ,fitter->GetPedestalRMS()) ; break ;
565         case 2: FillRawsData(kHGpedRMSMod2,cellX,cellZ,fitter->GetPedestalRMS()) ; break ;
566         case 3: FillRawsData(kHGpedRMSMod3,cellX,cellZ,fitter->GetPedestalRMS()) ; break ;
567         case 4: FillRawsData(kHGpedRMSMod4,cellX,cellZ,fitter->GetPedestalRMS()) ; break ;
568         }               
569         //if quality was evaluated, fill histo
570         if(strcmp(GetRecoParam()->EMCFitterVersion(),"v1")==0){
571           switch (module){
572           case 0: FillRawsData(kHGqualMod0,cellX,cellZ,fitter->GetSignalQuality()) ; break ;
573           case 1: FillRawsData(kHGqualMod1,cellX,cellZ,fitter->GetSignalQuality()) ; break ;
574           case 2: FillRawsData(kHGqualMod2,cellX,cellZ,fitter->GetSignalQuality()) ; break ;
575           case 3: FillRawsData(kHGqualMod3,cellX,cellZ,fitter->GetSignalQuality()) ; break ;
576           case 4: FillRawsData(kHGqualMod4,cellX,cellZ,fitter->GetSignalQuality()) ; break ;
577           }       
578         }
579         FillRawsData(kNmodHG,module) ; 
580         FillRawsData(kNRCUHG,RCUnum) ; 
581         FillRawsData(kHGtime,time) ;  
582         FillRawsData(kSpecHG,energy) ;
583         hgEtot+=energy ; 
584         hgNtot++ ;  
585         if(RCUnum==0)hgNtotM1RCU0++;
586         if(RCUnum==1)hgNtotM1RCU1++;
587         if(RCUnum==2)hgNtotM1RCU2++;
588         if(RCUnum==3)hgNtotM1RCU3++;
589         if(RCUnum==4)hgNtotM2RCU0++;
590         if(RCUnum==5)hgNtotM2RCU1++;
591         if(RCUnum==6)hgNtotM2RCU2++;
592         if(RCUnum==7)hgNtotM2RCU3++;
593         if(RCUnum==8)hgNtotM3RCU0++;
594         if(RCUnum==9)hgNtotM3RCU1++;
595         if(RCUnum==10)hgNtotM3RCU2++;
596         if(RCUnum==11)hgNtotM3RCU3++;
597       }
598     }  // End of NextChannel
599   } // End of NextDDL
600   delete fitter;
601
602   FillRawsData(kEtotLG,lgEtot) ; 
603   TParameter<double> * p;
604   p = dynamic_cast<TParameter<double>*>(GetParameterList()->
605                                         FindObject(Form("%s_%s_%s", GetName(), 
606                                                         AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), 
607                                                         GetRawsData(kEtotLG)->GetName()))) ; 
608   if (p) p->SetVal(lgEtot) ; 
609   FillRawsData(kEtotHG,hgEtot) ;  
610   p = dynamic_cast<TParameter<double>*>(GetParameterList()->
611                                         FindObject(Form("%s_%s_%s", GetName(), 
612                                                         AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), 
613                                                         GetRawsData(kEtotHG)->GetName()))) ; 
614   if (p) p->SetVal(hgEtot) ; 
615   FillRawsData(kNtotLG,lgNtot) ;
616   p = dynamic_cast<TParameter<double>*>(GetParameterList()->
617                                         FindObject(Form("%s_%s_%s", GetName(), 
618                                                         AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), 
619                                                         GetRawsData(kNtotLG)->GetName()))) ; 
620   if (p) p->SetVal(lgNtot) ; 
621   FillRawsData(kNtotHG,hgNtot) ;
622   p = dynamic_cast<TParameter<double>*>(GetParameterList()->
623                                         FindObject(Form("%s_%s_%s", 
624                                                         GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), 
625                                                         GetRawsData(kNtotHG)->GetName()))) ; 
626   if (p) p->SetVal(hgNtot) ; 
627
628   FillRawsData(kNtotHG_1_0,hgNtotM1RCU0) ;
629   p = dynamic_cast<TParameter<double>*>(GetParameterList()->
630                                         FindObject(Form("%s_%s_%s",
631                                                         GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(),
632                                                         GetRawsData(kNtotHG_1_0)->GetName()))) ;
633   if (p) p->SetVal(hgNtotM1RCU0) ;
634
635   FillRawsData(kNtotHG_1_1,hgNtotM1RCU1) ;
636   p = dynamic_cast<TParameter<double>*>(GetParameterList()->
637                                         FindObject(Form("%s_%s_%s",
638                                                         GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(),
639                                                         GetRawsData(kNtotHG_1_1)->GetName()))) ;
640   if (p) p->SetVal(hgNtotM1RCU1) ;
641
642   FillRawsData(kNtotHG_1_2,hgNtotM1RCU2) ;
643   p = dynamic_cast<TParameter<double>*>(GetParameterList()->
644                                         FindObject(Form("%s_%s_%s",
645                                                         GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(),
646                                                         GetRawsData(kNtotHG_1_2)->GetName()))) ;
647   if (p) p->SetVal(hgNtotM1RCU2) ;
648
649   FillRawsData(kNtotHG_1_3,hgNtotM1RCU3) ;
650   p = dynamic_cast<TParameter<double>*>(GetParameterList()->
651                                         FindObject(Form("%s_%s_%s",
652                                                         GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(),
653                                                         GetRawsData(kNtotHG_1_3)->GetName()))) ;
654   if (p) p->SetVal(hgNtotM1RCU3) ;
655
656   FillRawsData(kNtotHG_2_0,hgNtotM2RCU0) ;
657   p = dynamic_cast<TParameter<double>*>(GetParameterList()->
658                                         FindObject(Form("%s_%s_%s",
659                                                         GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(),
660                                                         GetRawsData(kNtotHG_2_0)->GetName()))) ;
661   if (p) p->SetVal(hgNtotM2RCU0) ;
662
663   FillRawsData(kNtotHG_2_1,hgNtotM2RCU1) ;
664   p = dynamic_cast<TParameter<double>*>(GetParameterList()->
665                                         FindObject(Form("%s_%s_%s",
666                                                         GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(),
667                                                         GetRawsData(kNtotHG_2_1)->GetName()))) ;
668   if (p) p->SetVal(hgNtotM2RCU1) ;
669
670   FillRawsData(kNtotHG_2_2,hgNtotM2RCU2) ;
671   p = dynamic_cast<TParameter<double>*>(GetParameterList()->
672                                         FindObject(Form("%s_%s_%s",
673                                                         GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(),
674                                                         GetRawsData(kNtotHG_2_2)->GetName()))) ;
675   if (p) p->SetVal(hgNtotM2RCU2) ;
676
677   FillRawsData(kNtotHG_2_3,hgNtotM2RCU3) ;
678   p = dynamic_cast<TParameter<double>*>(GetParameterList()->
679                                         FindObject(Form("%s_%s_%s",
680                                                         GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(),
681                                                         GetRawsData(kNtotHG_2_3)->GetName()))) ;
682   if (p) p->SetVal(hgNtotM2RCU3) ;
683
684   FillRawsData(kNtotHG_3_0,hgNtotM3RCU0) ;
685   p = dynamic_cast<TParameter<double>*>(GetParameterList()->
686                                         FindObject(Form("%s_%s_%s",
687                                                         GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(),
688                                                         GetRawsData(kNtotHG_3_0)->GetName()))) ;
689   if (p) p->SetVal(hgNtotM3RCU0) ;
690
691   FillRawsData(kNtotHG_3_1,hgNtotM3RCU1) ;
692   p = dynamic_cast<TParameter<double>*>(GetParameterList()->
693                                         FindObject(Form("%s_%s_%s",
694                                                         GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(),
695                                                         GetRawsData(kNtotHG_3_1)->GetName()))) ;
696   if (p) p->SetVal(hgNtotM3RCU1) ;
697
698   FillRawsData(kNtotHG_3_2,hgNtotM3RCU2) ;
699   p = dynamic_cast<TParameter<double>*>(GetParameterList()->
700                                         FindObject(Form("%s_%s_%s",
701                                                         GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(),
702                                                         GetRawsData(kNtotHG_3_2)->GetName()))) ;
703   if (p) p->SetVal(hgNtotM3RCU2) ;
704
705   FillRawsData(kNtotHG_3_3,hgNtotM3RCU3) ;
706   p = dynamic_cast<TParameter<double>*>(GetParameterList()->
707                                         FindObject(Form("%s_%s_%s",
708                                                         GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(),
709                                                         GetRawsData(kNtotHG_3_3)->GetName()))) ;
710   if (p) p->SetVal(hgNtotM3RCU3) ;
711
712   //
713   IncEvCountCycleRaws();
714   IncEvCountTotalRaws();
715   //
716 }
717
718 //____________________________________________________________________________
719 void AliPHOSQADataMakerRec::MakeDigits()
720 {
721   // makes data from Digits
722   
723   if ( ! GetDigitsData(kDigitsMul) ) InitDigits() ;
724   FillDigitsData(kDigitsMul,fDigitsArray->GetEntriesFast()) ; 
725   TIter next(fDigitsArray) ; 
726   AliPHOSDigit * digit ; 
727   while ( (digit = dynamic_cast<AliPHOSDigit *>(next())) ) {
728     FillDigitsData(kDigits, digit->GetEnergy()) ;
729   }  
730   //
731 }
732
733 //____________________________________________________________________________
734 void AliPHOSQADataMakerRec::MakeDigits(TTree * digitTree)
735 {
736   // makes data from Digit Tree
737   if (fDigitsArray) 
738     fDigitsArray->Clear() ; 
739   else 
740     fDigitsArray = new TClonesArray("AliPHOSDigit", 1000) ; 
741   
742   TBranch * branch = digitTree->GetBranch("PHOS") ;
743   if ( ! branch ) {AliWarning("PHOS branch in Digit Tree not found"); return;} 
744   branch->SetAddress(&fDigitsArray) ;
745   branch->GetEntry(0) ; 
746   MakeDigits() ; 
747   //
748   IncEvCountCycleDigits();
749   IncEvCountTotalDigits();
750   //
751 }
752
753 //____________________________________________________________________________
754 void AliPHOSQADataMakerRec::MakeRecPoints(TTree * clustersTree)
755 {
756   {
757     // makes data from RecPoints
758     TBranch *emcbranch = clustersTree->GetBranch("PHOSEmcRP");
759     if (!emcbranch) { 
760       AliError("can't get the branch with the PHOS EMC clusters !");
761       return;
762     }
763     
764     TObjArray * emcrecpoints = new TObjArray(100) ;
765     emcbranch->SetAddress(&emcrecpoints);
766     emcbranch->GetEntry(0);
767     
768     FillRecPointsData(kRPNtot,emcrecpoints->GetEntriesFast()) ; 
769     TIter next(emcrecpoints) ; 
770     AliPHOSEmcRecPoint * rp ; 
771     Double_t eTot = 0. ; 
772     while ( (rp = static_cast<AliPHOSEmcRecPoint *>(next())) ) {
773       FillRecPointsData(kRPSpec, rp->GetEnergy()) ;
774       Int_t mod = rp->GetPHOSMod() ;
775       TVector3 pos ;
776       rp->GetLocalPosition(pos) ;
777       switch(mod){
778         case 1: FillRecPointsData(kRPmod1,pos.X(),pos.Z()) ; break ;
779         case 2: FillRecPointsData(kRPmod2,pos.X(),pos.Z()) ; break ;
780         case 3: FillRecPointsData(kRPmod3,pos.X(),pos.Z()) ; break ;
781         case 4: FillRecPointsData(kRPmod4,pos.X(),pos.Z()) ; break ;
782         case 5: FillRecPointsData(kRPmod5,pos.X(),pos.Z()) ; break ;
783       }
784       eTot+= rp->GetEnergy() ;
785     }
786     FillRecPointsData(kRPEtot,eTot) ;
787     emcrecpoints->Delete();
788     delete emcrecpoints;
789   }
790   {
791     TBranch *cpvbranch = clustersTree->GetBranch("PHOSCpvRP");
792     if (!cpvbranch) { 
793       AliError("can't get the branch with the PHOS CPV clusters !");
794       return;
795     }
796     TObjArray *cpvrecpoints = new TObjArray(100) ;
797     cpvbranch->SetAddress(&cpvrecpoints);
798     cpvbranch->GetEntry(0);
799     
800     FillRecPointsData(kRPNcpv,cpvrecpoints->GetEntriesFast()) ; 
801     cpvrecpoints->Delete();
802     delete cpvrecpoints;
803   }
804   //
805   IncEvCountCycleRecPoints();
806   IncEvCountTotalRecPoints();
807   //
808 }
809
810 //____________________________________________________________________________ 
811 void AliPHOSQADataMakerRec::StartOfDetectorCycle()
812 {
813   //Detector specific actions at start of cycle
814   
815 }