1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
20 Produces the data needed to calculate the quality assurance.
21 All data must be mergeable objects.
22 Y. Schutz CERN July 2007
25 // --- ROOT system ---
26 #include <TClonesArray.h>
31 #include <TParameter.h>
33 // --- Standard library ---
35 // --- AliRoot header files ---
36 #include "AliCaloRawStreamV3.h"
37 #include "AliESDCaloCluster.h"
38 #include "AliESDEvent.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"
53 ClassImp(AliPHOSQADataMakerRec)
55 //____________________________________________________________________________
56 AliPHOSQADataMakerRec::AliPHOSQADataMakerRec() :
57 AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kPHOS), "PHOS Quality Assurance Data Maker")
62 //____________________________________________________________________________
63 AliPHOSQADataMakerRec::AliPHOSQADataMakerRec(const AliPHOSQADataMakerRec& qadm) :
67 SetName((const char*)qadm.GetName()) ;
68 SetTitle((const char*)qadm.GetTitle());
71 //__________________________________________________________________
72 AliPHOSQADataMakerRec& AliPHOSQADataMakerRec::operator = (const AliPHOSQADataMakerRec& qadm )
75 this->~AliPHOSQADataMakerRec();
76 new(this) AliPHOSQADataMakerRec(qadm);
80 //____________________________________________________________________________
81 void AliPHOSQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
83 //Detector specific actions at end of cycle
84 ResetEventTrigClasses();
86 for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
87 if (! IsValidEventSpecie(specie, list)) continue;
88 SetEventSpecie(AliRecoParam::ConvertIndex(specie));
90 for (int itc=-1;itc<GetNTrigClasses();itc++) { // RS: loop over all trigger clones
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 } // RS: loop over all trigger clones
104 // do the QA checking
105 AliQAChecker::Instance()->Run(AliQAv1::kPHOS, task, list) ;
108 //____________________________________________________________________________
109 void AliPHOSQADataMakerRec::InitESDs()
111 //Create histograms to controll ESD
113 const Bool_t expert = kTRUE ;
114 const Bool_t image = kTRUE ;
116 TH1F * h1 = new TH1F("hESDPhosSpectrum", "ESDs spectrum in PHOS; Energy [MeV];Counts" , 200, 0., 20.) ;
118 Add2ESDsList(h1, kESDSpec, !expert, image) ;
120 TH1I * h2 = new TH1I("hESDPhosMul", "ESDs multiplicity distribution in PHOS; # of clusters;Counts", 100, 0, 100 ) ;
122 Add2ESDsList(h2, kESDNtot, !expert, image) ;
124 TH1F * h3 = new TH1F("hESDPhosEtot", "ESDs total energy;Energy [MeV];Counts" , 2000, 0, 200.) ;
126 Add2ESDsList(h3, kESDEtot, !expert, image) ; //Expert histo
128 TH1F * h4 = new TH1F("hESDpid", "ESDs PID distribution in PHOS;Particle Id;Counts" , 100, 0., 1.) ;
130 Add2ESDsList(h4, kESDpid, !expert, image) ; //Expert histo
132 ClonePerTrigClass(AliQAv1::kESDS); // this should be the last line
135 //____________________________________________________________________________
136 void AliPHOSQADataMakerRec::InitDigits()
138 // create Digits histograms in Digits subdir
139 const Bool_t expert = kTRUE ;
140 const Bool_t image = kTRUE ;
141 TH1I * h0 = new TH1I("hPhosDigits", "Digits amplitude distribution in PHOS;Amplitude [ADC counts];Counts", 500, 0, 1000) ;
143 Add2DigitsList(h0, kDigits, !expert, image) ;
144 TH1I * h1 = new TH1I("hPhosDigitsMul", "Digits multiplicity distribution in PHOS;# of Digits;Entries", 2000, 0, 10000) ;
146 Add2DigitsList(h1, kDigitsMul, !expert, image) ;
148 ClonePerTrigClass(AliQAv1::kDIGITS); // this should be the last line
151 //____________________________________________________________________________
152 void AliPHOSQADataMakerRec::InitRecPoints()
154 // create Reconstructed Points histograms in RecPoints subdir
155 const Bool_t expert = kTRUE ;
156 const Bool_t image = kTRUE ;
158 TH2I * h0 = new TH2I("hRpPHOSxyMod1","RecPoints Rows x Columns for PHOS module 1;Row #;Column #", 64, -72., 72., 56, -63., 63.) ;
159 Add2RecPointsList(h0,kRPmod1, expert, !image) ;
160 TH2I * h1 = new TH2I("hRpPHOSxyMod2","RecPoints Rows x Columns for PHOS module 2;Row #;Column #", 64, -72., 72., 56, -63., 63.) ;
161 Add2RecPointsList(h1,kRPmod2, expert, !image) ;
162 TH2I * h2 = new TH2I("hRpPHOSxyMod3","RecPoints Rows x Columns for PHOS module 3;Row #;Column #", 64, -72., 72., 56, -63., 63.) ;
163 Add2RecPointsList(h2,kRPmod3, expert, !image) ;
164 TH2I * h3 = new TH2I("hRpPHOSxyMod4","RecPoints Rows x Columns for PHOS module 4;Row #;Column #", 64, -72., 72., 56, -63., 63.) ;
165 Add2RecPointsList(h3,kRPmod4, expert, !image) ;
166 TH2I * h4 = new TH2I("hRpPHOSxyMod5","RecPoints Rows x Columns for PHOS module 5;Row #;Column #", 64, -72., 72., 56, -63., 63.) ;
167 Add2RecPointsList(h4,kRPmod5, expert, !image) ;
169 TH1F * h5 = new TH1F("hEmcPhosRecPointsSpectrum", "EMC RecPoints spectrum in PHOS;Energy [MeV];Counts", 2000, 0., 20.) ;
171 Add2RecPointsList(h5, kRPSpec, !expert, image) ;
173 TH1I * h6 = new TH1I("hEmcPhosRecPointsMul", "EMC RecPoints multiplicity distribution in PHOS;# of EMC Clusters;Entries", 100, 0, 100) ;
175 Add2RecPointsList(h6, kRPNtot, !expert, image) ;
177 TH1I * h7 = new TH1I("hEmcPhosRecPointsEtot", "EMC RecPoints Etot;Energy [MeV];Counts", 200, 0, 200.) ;
179 Add2RecPointsList(h7, kRPEtot, !expert, image) ;
181 TH1I * h8 = new TH1I("hCpvPhosRecPointsMul", "CPV RecPoints multiplicity distribution in PHOS;# of CPV clusters;Counts", 100, 0, 100) ;
183 Add2RecPointsList(h8, kRPNcpv, !expert, image) ;
185 ClonePerTrigClass(AliQAv1::kRECPOINTS); // this should be the last line
188 //____________________________________________________________________________
189 void AliPHOSQADataMakerRec::InitRaws()
191 // create Raws histograms in Raws subdir
192 const Bool_t expert = kTRUE ;
193 const Bool_t saveCorr = kTRUE ;
194 const Bool_t image = kTRUE ;
196 TH2I * h0 = new TH2I("hHighPHOSxyMod0","High Gain in PHOS module 0", 64, 0, 64, 56, 0, 56) ;
197 h0->SetXTitle("x, cells"); h0->SetYTitle("z, cells");
198 Add2RawsList(h0,kHGmod0, expert, !image, !saveCorr) ;
199 TH2I * h1 = new TH2I("hHighPHOSxyMod1","High Gain in PHOS module 1", 64, 0, 64, 56, 0, 56) ;
200 h1->SetXTitle("x, cells"); h1->SetYTitle("z, cells");
201 Add2RawsList(h1,kHGmod1, expert, !image, !saveCorr) ;
202 TH2I * h2 = new TH2I("hHighPHOSxyMod2","High Gain in PHOS module 2", 64, 0, 64, 56, 0, 56) ;
203 h2->SetXTitle("x, cells"); h2->SetYTitle("z, cells");
204 Add2RawsList(h2,kHGmod2, expert, !image, !saveCorr) ;
205 TH2I * h3 = new TH2I("hHighPHOSxyMod3","High Gain in PHOS module 3", 64, 0, 64, 56, 0, 56) ;
206 h3->SetXTitle("x, cells"); h3->SetYTitle("z, cells");
207 Add2RawsList(h3,kHGmod3, expert, !image, !saveCorr) ;
208 TH2I * h4 = new TH2I("hHighPHOSxyMod4","High Gain in PHOS module 4", 64, 0, 64, 56, 0, 56) ;
209 h4->SetXTitle("x, cells"); h4->SetYTitle("z, cells");
210 Add2RawsList(h4,kHGmod4, expert, !image, !saveCorr) ;
212 TH2I * h5 = new TH2I("hLowPHOSxyMod0","Low Gain in PHOS module 0", 64, 0, 64, 56, 0, 56) ;
213 h5->SetXTitle("x, cells"); h5->SetYTitle("z, cells");
214 Add2RawsList(h5,kLGmod0, expert, !image, !saveCorr) ;
215 TH2I * h6 = new TH2I("hLowPHOSxyMod1","Low Gain in PHOS module 1", 64, 0, 64, 56, 0, 56) ;
216 h6->SetXTitle("x, cells"); h6->SetYTitle("z, cells");
217 Add2RawsList(h6,kLGmod1, expert, !image, !saveCorr) ;
218 TH2I * h7 = new TH2I("hLowPHOSxyMod2","Low Gain in PHOS module 2", 64, 0, 64, 56, 0, 56) ;
219 h7->SetXTitle("x, cells"); h7->SetYTitle("z, cells");
220 Add2RawsList(h7,kLGmod2, expert, !image, !saveCorr) ;
221 TH2I * h8 = new TH2I("hLowPHOSxyMod3","Low Gain in PHOS module 3", 64, 0, 64, 56, 0, 56) ;
222 h8->SetXTitle("x, cells"); h8->SetYTitle("z, cells");
223 Add2RawsList(h8,kLGmod3, expert, !image, !saveCorr) ;
224 TH2I * h9 = new TH2I("hLowPHOSxyMod4","Low Gain in PHOS module 4", 64, 0, 64, 56, 0, 56) ;
225 h9->SetXTitle("x, cells"); h9->SetYTitle("z, cells");
226 Add2RawsList(h9,kLGmod4, expert, !image, !saveCorr) ;
228 TH1I * h10 = new TH1I("hLowPhosModules", "Low Gain Hits in EMCA PHOS modules", 5, 0, 5) ;
229 h10->SetXTitle("Module number");
230 Add2RawsList(h10, kNmodLG, !expert, image, !saveCorr) ;
231 TH1I * h11 = new TH1I("hHighPhosModules", "High Gain Hits in EMCA PHOS modules",5, 0, 5) ;
232 h11->SetXTitle("Module number");
233 Add2RawsList(h11, kNmodHG, !expert, image, !saveCorr) ;
235 TH1F * h12 = new TH1F("hLowPhosRawtime" , "Low Gain Time of raw hits in PHOS" , 500, -50., 200.) ;
236 h12->SetXTitle("Time [samples]");
238 Add2RawsList(h12, kLGtime, expert, !image, !saveCorr) ;
239 TH1F * h13 = new TH1F("hHighPhosRawtime", "High Gain Time of raw hits in PHOS", 500, -50., 200.) ;
240 h13->SetXTitle("Time [samples]");
242 Add2RawsList(h13, kHGtime, expert, !image, !saveCorr) ;
244 TH1F * h14 = new TH1F("hLowPhosRawEnergy" , "Low Gain Energy of raw hits in PHOS" , 512, 0., 1024.) ;
245 h14->SetXTitle("Energy [ADC counts]");
247 Add2RawsList(h14, kSpecLG, !expert, image, !saveCorr) ;
248 TH1F * h15 = new TH1F("hHighPhosRawEnergy", "High Gain Energy of raw hits in PHOS", 512, 0., 1024.) ;
249 h15->SetXTitle("Energy [ADC counts]");
251 Add2RawsList(h15, kSpecHG, !expert, image, !saveCorr) ;
253 TH1F * h16 = new TH1F("hLowNtot" , "Low Gain Total Number of raw hits in PHOS" , 500, 0., 5000.) ;
254 h16->SetXTitle("Number of hits");
256 Add2RawsList(h16, kNtotLG, !expert, image, saveCorr) ;
257 TH1F * h17 = new TH1F("hHighNtot", "High Gain Total Number of raw hits in PHOS", 500, 0., 5000.) ;
258 h17->SetXTitle("Number of hits");
260 Add2RawsList(h17, kNtotHG, !expert, image, saveCorr) ;
262 TH1F * h18 = new TH1F("hLowEtot" , "Low Gain Total Energy of raw hits in PHOS" , 500, 0., 100000.) ;
263 h18->SetXTitle("Energy [ADC counts]");
265 Add2RawsList(h18, kEtotLG, !expert, image, saveCorr) ;
266 TH1F * h19 = new TH1F("hHighEtot", "High Gain Total Energy of raw hits in PHOS", 500, 0., 100000.) ;
267 h19->SetXTitle("Energy [ADC counts]");
269 Add2RawsList(h19, kEtotHG, !expert, image, saveCorr) ;
271 TH2F * h20 = new TH2F("hQualHGxyMod0","High Gain signal quality in module 0", 64, 0, 64, 56, 0, 56) ;
272 h20->SetXTitle("x, cells"); h20->SetYTitle("z, cells");
273 Add2RawsList(h20,kHGqualMod0, expert, !image, !saveCorr) ;
274 TH2F * h21 = new TH2F("hQualHGxyMod1","High Gain signal quality in module 1", 64, 0, 64, 56, 0, 56) ;
275 h21->SetXTitle("x, cells"); h21->SetYTitle("z, cells");
276 Add2RawsList(h21,kHGqualMod1, expert, !image, !saveCorr) ;
277 TH2F * h22 = new TH2F("hQualHGxyMod2","High Gain signal quality in module 2", 64, 0, 64, 56, 0, 56) ;
278 h22->SetXTitle("x, cells"); h22->SetYTitle("z, cells");
279 Add2RawsList(h22,kHGqualMod2, expert, !image, !saveCorr) ;
280 TH2F * h23 = new TH2F("hQualHGxyMod3","High Gain signal quality in module 3", 64, 0, 64, 56, 0, 56) ;
281 h23->SetXTitle("x, cells"); h23->SetYTitle("z, cells");
282 Add2RawsList(h23,kHGqualMod3, expert, !image, !saveCorr) ;
283 TH2F * h24 = new TH2F("hQualHGxyMod4","High Gain signal quality in module 4", 64, 0, 64, 56, 0, 56) ;
284 h24->SetXTitle("x, cells"); h24->SetYTitle("z, cells");
285 Add2RawsList(h24,kHGqualMod4, expert, !image, !saveCorr) ;
287 TH2F * h25 = new TH2F("hQualLGxyMod0","Low Gain signal quality in module 0", 64, 0, 64, 56, 0, 56) ;
288 h25->SetXTitle("x, cells"); h25->SetYTitle("z, cells");
289 Add2RawsList(h25,kLGqualMod0, expert, !image, !saveCorr) ;
290 TH2F * h26 = new TH2F("hQualLGxyMod1","Low Gain signal quality in module 1", 64, 0, 64, 56, 0, 56) ;
291 h26->SetXTitle("x, cells"); h26->SetYTitle("z, cells");
292 Add2RawsList(h26,kLGqualMod1, expert, !image, !saveCorr) ;
293 TH2F * h27 = new TH2F("hQualLGxyMod2","Low Gain signal quality in module 2", 64, 0, 64, 56, 0, 56) ;
294 h27->SetXTitle("x, cells"); h27->SetYTitle("z, cells");
295 Add2RawsList(h27,kLGqualMod2, expert, !image, !saveCorr) ;
296 TH2F * h28 = new TH2F("hQualLGxyMod3","Low Gain signal quality in module 3", 64, 0, 64, 56, 0, 56) ;
297 h28->SetXTitle("x, cells"); h28->SetYTitle("z, cells");
298 Add2RawsList(h28,kLGqualMod3, expert, !image, !saveCorr) ;
299 TH2F * h29 = new TH2F("hQualLGxyMod4","Low Gain signal quality in module 4", 64, 0, 64, 56, 0, 56) ;
300 h29->SetXTitle("x, cells"); h29->SetYTitle("z, cells");
301 Add2RawsList(h29,kLGqualMod4, expert, !image, !saveCorr) ;
303 TH1F * h30 = new TH1F("hLGpedRMS" ,"Low Gain pedestal RMS" ,200,0.,20.) ;
304 h30->SetXTitle("RMS [ADC counts]");
306 Add2RawsList(h30,kLGpedRMS, expert, !image, !saveCorr) ;
307 TH1F * h31 = new TH1F("hHGpedRMS" ,"High Gain pedestal RMS",200,0.,20.) ;
308 h31->SetXTitle("RMS [ADC counts]");
310 Add2RawsList(h31,kHGpedRMS, expert, !image, !saveCorr) ;
312 TH2F * h32 = new TH2F("hpedRMSHGxyMod0","High Gain pedestal RMS in module 0", 64, 0, 64, 56, 0, 56) ;
313 h32->SetXTitle("x, cells"); h32->SetYTitle("z, cells");
314 Add2RawsList(h32,kHGpedRMSMod0, expert, !image, !saveCorr) ;
315 TH2F * h33 = new TH2F("hpedRMSHGxyMod1","High Gain pedestal RMS in module 1", 64, 0, 64, 56, 0, 56) ;
316 h33->SetXTitle("x, cells"); h33->SetYTitle("z, cells");
317 Add2RawsList(h33,kHGpedRMSMod1, expert, !image, !saveCorr) ;
318 TH2F * h34 = new TH2F("hpedRMSHGxyMod2","High Gain pedestal RMS in module 2", 64, 0, 64, 56, 0, 56) ;
319 h34->SetXTitle("x, cells"); h34->SetYTitle("z, cells");
320 Add2RawsList(h34,kHGpedRMSMod2, expert, !image, !saveCorr) ;
321 TH2F * h35 = new TH2F("hpedRMSHGxyMod3","High Gain pedestal RMS in module 3", 64, 0, 64, 56, 0, 56) ;
322 h35->SetXTitle("x, cells"); h35->SetYTitle("z, cells");
323 Add2RawsList(h35,kHGpedRMSMod3, expert, !image, !saveCorr) ;
324 TH2F * h36 = new TH2F("hpedRMSHGxyMod4","High Gain pedestal RMS in module 4", 64, 0, 64, 56, 0, 56) ;
325 h36->SetXTitle("x, cells"); h36->SetYTitle("z, cells");
326 Add2RawsList(h36,kHGpedRMSMod4, expert, !image, !saveCorr) ;
328 TH2F * h37 = new TH2F("hpedRMSLGxyMod0","Low Gain pedestal RMS in module 0", 64, 0, 64, 56, 0, 56) ;
329 h37->SetXTitle("x, cells"); h37->SetYTitle("z, cells");
330 Add2RawsList(h37,kLGpedRMSMod0, expert, !image, !saveCorr) ;
331 TH2F * h38 = new TH2F("hpedRMSLGxyMod1","Low Gain pedestal RMS in module 1", 64, 0, 64, 56, 0, 56) ;
332 h38->SetXTitle("x, cells"); h38->SetYTitle("z, cells");
333 Add2RawsList(h38,kLGpedRMSMod1, expert, !image, !saveCorr) ;
334 TH2F * h39 = new TH2F("hpedRMSLGxyMod2","Low Gain pedestal RMS in module 2", 64, 0, 64, 56, 0, 56) ;
335 h39->SetXTitle("x, cells"); h39->SetYTitle("z, cells");
336 Add2RawsList(h39,kLGpedRMSMod2, expert, !image, !saveCorr) ;
337 TH2F * h40 = new TH2F("hpedRMSLGxyMod3","Low Gain pedestal RMS in module 3", 64, 0, 64, 56, 0, 56) ;
338 h40->SetXTitle("x, cells"); h40->SetYTitle("z, cells");
339 Add2RawsList(h40,kLGpedRMSMod3, expert, !image, !saveCorr) ;
340 TH2F * h41 = new TH2F("hpedRMSLGxyMod4","Low Gain pedestal RMS in module 4", 64, 0, 64, 56, 0, 56) ;
341 h41->SetXTitle("x, cells"); h41->SetYTitle("z, cells");
342 Add2RawsList(h41,kLGpedRMSMod4, expert, !image, !saveCorr) ;
344 ClonePerTrigClass(AliQAv1::kRAWS); // this should be the last line
347 //____________________________________________________________________________
348 void AliPHOSQADataMakerRec::MakeESDs(AliESDEvent * esd)
350 // make QA data from ESDs
353 for ( Int_t index = 0; index < esd->GetNumberOfCaloClusters() ; index++ ) {
354 AliESDCaloCluster * clu = esd->GetCaloCluster(index) ;
355 if( clu->IsPHOS() ) {
356 FillESDsData(kESDSpec,clu->E()) ;
357 const Double_t * pid = clu->GetPID() ;
358 FillESDsData(kESDpid,pid[AliPID::kPhoton]) ;
363 FillESDsData(kESDNtot,nTot) ;
364 FillESDsData(kESDEtot,eTot) ;
366 IncEvCountCycleESDs();
367 IncEvCountTotalESDs();
371 //____________________________________________________________________________
372 void AliPHOSQADataMakerRec::MakeRaws(AliRawReader* rawReader)
374 //Fill prepared histograms with Raw digit properties
378 const TObjArray* maps = AliPHOSRecoParam::GetMappings();
379 if(!maps) AliFatal("Cannot retrieve ALTRO mappings!!");
381 AliAltroMapping *mapping[20];
382 for(Int_t i = 0; i < 20; i++) {
383 mapping[i] = (AliAltroMapping*)maps->At(i);
386 AliCaloRawStreamV3 fRawStream(rawReader,"PHOS",mapping);
388 AliPHOSRawFitterv0 * fitter ;
389 if (strcmp(GetRecoParam()->EMCFitterVersion(),"v1")==0)
390 fitter=new AliPHOSRawFitterv1();
391 else if(strcmp(GetRecoParam()->EMCFitterVersion(),"v2")==0)
392 fitter=new AliPHOSRawFitterv2();
394 fitter=new AliPHOSRawFitterv0();
400 while (fRawStream.NextDDL()) {
401 while (fRawStream.NextChannel()) {
402 Int_t module = fRawStream.GetModule();
403 Int_t cellX = fRawStream.GetCellX();
404 Int_t cellZ = fRawStream.GetCellZ();
405 Int_t caloFlag = fRawStream.GetCaloFlag(); // 0=LG, 1=HG, 2=TRU
407 if(caloFlag!=0 && caloFlag!=1) continue; //TRU data!
409 fitter->SetChannelGeo(module+1,cellX+1,cellZ+1,caloFlag);
411 if(fitter->GetAmpOffset()==0 && fitter->GetAmpThreshold()==0) {
412 Short_t altroCFG1 = fRawStream.GetAltroCFG1();
413 Bool_t ZeroSuppressionEnabled = (altroCFG1 >> 15) & 0x1;
414 if(ZeroSuppressionEnabled) {
415 Short_t offset = (altroCFG1 >> 10) & 0xf;
416 Short_t threshold = altroCFG1 & 0x3ff;
417 fitter->SubtractPedestals(kFALSE);
418 fitter->SetAmpOffset(offset);
419 fitter->SetAmpThreshold(threshold);
422 fitter->SubtractPedestals(kTRUE);
426 while (fRawStream.NextBunch()) {
428 if (nBunches > 1) continue;
429 const UShort_t *sig = fRawStream.GetSignals();
430 Int_t sigStart = fRawStream.GetStartTimeBin();
431 Int_t sigLength = fRawStream.GetBunchLength();
432 fitter->Eval(sig,sigStart,sigLength);
433 } // End of NextBunch()
435 Double_t energy = fitter->GetEnergy() ;
436 Double_t time = fitter->GetTime() ;
438 if (caloFlag == 0) { // LG
439 FillRawsData(kLGpedRMS,fitter->GetPedestalRMS()) ;
441 case 0: FillRawsData(kLGmod0,cellX,cellZ) ; break ;
442 case 1: FillRawsData(kLGmod1,cellX,cellZ) ; break ;
443 case 2: FillRawsData(kLGmod2,cellX,cellZ) ; break ;
444 case 3: FillRawsData(kLGmod3,cellX,cellZ) ; break ;
445 case 4: FillRawsData(kLGmod4,cellX,cellZ) ; break ;
448 case 0: FillRawsData(kLGpedRMSMod0,cellX,cellZ,fitter->GetPedestalRMS()) ; break ;
449 case 1: FillRawsData(kLGpedRMSMod1,cellX,cellZ,fitter->GetPedestalRMS()) ; break ;
450 case 2: FillRawsData(kLGpedRMSMod2,cellX,cellZ,fitter->GetPedestalRMS()) ; break ;
451 case 3: FillRawsData(kLGpedRMSMod3,cellX,cellZ,fitter->GetPedestalRMS()) ; break ;
452 case 4: FillRawsData(kLGpedRMSMod4,cellX,cellZ,fitter->GetPedestalRMS()) ; break ;
454 //if quality was evaluated, fill histo
455 if(strcmp(GetRecoParam()->EMCFitterVersion(),"v1")==0){
457 case 0: FillRawsData(kLGqualMod0,cellX,cellZ,fitter->GetSignalQuality()) ; break ;
458 case 1: FillRawsData(kLGqualMod1,cellX,cellZ,fitter->GetSignalQuality()) ; break ;
459 case 2: FillRawsData(kLGqualMod2,cellX,cellZ,fitter->GetSignalQuality()) ; break ;
460 case 3: FillRawsData(kLGqualMod3,cellX,cellZ,fitter->GetSignalQuality()) ; break ;
461 case 4: FillRawsData(kLGqualMod4,cellX,cellZ,fitter->GetSignalQuality()) ; break ;
464 FillRawsData(kNmodLG,module) ;
465 FillRawsData(kLGtime,time) ;
466 FillRawsData(kSpecLG,energy) ;
470 else if (caloFlag == 1) { // HG
471 FillRawsData(kHGpedRMS,fitter->GetPedestalRMS()) ;
473 case 0: FillRawsData(kHGmod0,cellX,cellZ) ; break ;
474 case 1: FillRawsData(kHGmod1,cellX,cellZ) ; break ;
475 case 2: FillRawsData(kHGmod2,cellX,cellZ) ; break ;
476 case 3: FillRawsData(kHGmod3,cellX,cellZ) ; break ;
477 case 4: FillRawsData(kHGmod4,cellX,cellZ) ; break ;
480 case 0: FillRawsData(kHGpedRMSMod0,cellX,cellZ,fitter->GetPedestalRMS()) ; break ;
481 case 1: FillRawsData(kHGpedRMSMod1,cellX,cellZ,fitter->GetPedestalRMS()) ; break ;
482 case 2: FillRawsData(kHGpedRMSMod2,cellX,cellZ,fitter->GetPedestalRMS()) ; break ;
483 case 3: FillRawsData(kHGpedRMSMod3,cellX,cellZ,fitter->GetPedestalRMS()) ; break ;
484 case 4: FillRawsData(kHGpedRMSMod4,cellX,cellZ,fitter->GetPedestalRMS()) ; break ;
486 //if quality was evaluated, fill histo
487 if(strcmp(GetRecoParam()->EMCFitterVersion(),"v1")==0){
489 case 0: FillRawsData(kHGqualMod0,cellX,cellZ,fitter->GetSignalQuality()) ; break ;
490 case 1: FillRawsData(kHGqualMod1,cellX,cellZ,fitter->GetSignalQuality()) ; break ;
491 case 2: FillRawsData(kHGqualMod2,cellX,cellZ,fitter->GetSignalQuality()) ; break ;
492 case 3: FillRawsData(kHGqualMod3,cellX,cellZ,fitter->GetSignalQuality()) ; break ;
493 case 4: FillRawsData(kHGqualMod4,cellX,cellZ,fitter->GetSignalQuality()) ; break ;
496 FillRawsData(kNmodHG,module) ;
497 FillRawsData(kHGtime,time) ;
498 FillRawsData(kSpecHG,energy) ;
502 } // End of NextChannel
506 FillRawsData(kEtotLG,lgEtot) ;
507 TParameter<double> * p;
508 p = dynamic_cast<TParameter<double>*>(GetParameterList()->
509 FindObject(Form("%s_%s_%s", GetName(),
510 AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(),
511 GetRawsData(kEtotLG)->GetName()))) ;
512 if (p) p->SetVal(lgEtot) ;
513 FillRawsData(kEtotHG,hgEtot) ;
514 p = dynamic_cast<TParameter<double>*>(GetParameterList()->
515 FindObject(Form("%s_%s_%s", GetName(),
516 AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(),
517 GetRawsData(kEtotHG)->GetName()))) ;
518 if (p) p->SetVal(hgEtot) ;
519 FillRawsData(kNtotLG,lgNtot) ;
520 p = dynamic_cast<TParameter<double>*>(GetParameterList()->
521 FindObject(Form("%s_%s_%s", GetName(),
522 AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(),
523 GetRawsData(kNtotLG)->GetName()))) ;
524 if (p) p->SetVal(lgNtot) ;
525 FillRawsData(kNtotHG,hgNtot) ;
526 p = dynamic_cast<TParameter<double>*>(GetParameterList()->
527 FindObject(Form("%s_%s_%s",
528 GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(),
529 GetRawsData(kNtotHG)->GetName()))) ;
530 if (p) p->SetVal(hgNtot) ;
532 IncEvCountCycleRaws();
533 IncEvCountTotalRaws();
537 //____________________________________________________________________________
538 void AliPHOSQADataMakerRec::MakeDigits()
540 // makes data from Digits
542 if ( ! GetDigitsData(kDigitsMul) ) InitDigits() ;
543 FillDigitsData(kDigitsMul,fDigitsArray->GetEntriesFast()) ;
544 TIter next(fDigitsArray) ;
545 AliPHOSDigit * digit ;
546 while ( (digit = dynamic_cast<AliPHOSDigit *>(next())) ) {
547 FillDigitsData(kDigits, digit->GetEnergy()) ;
552 //____________________________________________________________________________
553 void AliPHOSQADataMakerRec::MakeDigits(TTree * digitTree)
555 // makes data from Digit Tree
557 fDigitsArray->Clear() ;
559 fDigitsArray = new TClonesArray("AliPHOSDigit", 1000) ;
561 TBranch * branch = digitTree->GetBranch("PHOS") ;
562 if ( ! branch ) {AliWarning("PHOS branch in Digit Tree not found"); return;}
563 branch->SetAddress(&fDigitsArray) ;
564 branch->GetEntry(0) ;
567 IncEvCountCycleDigits();
568 IncEvCountTotalDigits();
572 //____________________________________________________________________________
573 void AliPHOSQADataMakerRec::MakeRecPoints(TTree * clustersTree)
576 // makes data from RecPoints
577 TBranch *emcbranch = clustersTree->GetBranch("PHOSEmcRP");
579 AliError("can't get the branch with the PHOS EMC clusters !");
583 TObjArray * emcrecpoints = new TObjArray(100) ;
584 emcbranch->SetAddress(&emcrecpoints);
585 emcbranch->GetEntry(0);
587 FillRecPointsData(kRPNtot,emcrecpoints->GetEntriesFast()) ;
588 TIter next(emcrecpoints) ;
589 AliPHOSEmcRecPoint * rp ;
591 while ( (rp = static_cast<AliPHOSEmcRecPoint *>(next())) ) {
592 FillRecPointsData(kRPSpec, rp->GetEnergy()) ;
593 Int_t mod = rp->GetPHOSMod() ;
595 rp->GetLocalPosition(pos) ;
597 case 1: FillRecPointsData(kRPmod1,pos.X(),pos.Z()) ; break ;
598 case 2: FillRecPointsData(kRPmod2,pos.X(),pos.Z()) ; break ;
599 case 3: FillRecPointsData(kRPmod3,pos.X(),pos.Z()) ; break ;
600 case 4: FillRecPointsData(kRPmod4,pos.X(),pos.Z()) ; break ;
601 case 5: FillRecPointsData(kRPmod5,pos.X(),pos.Z()) ; break ;
603 eTot+= rp->GetEnergy() ;
605 FillRecPointsData(kRPEtot,eTot) ;
606 emcrecpoints->Delete();
610 TBranch *cpvbranch = clustersTree->GetBranch("PHOSCpvRP");
612 AliError("can't get the branch with the PHOS CPV clusters !");
615 TObjArray *cpvrecpoints = new TObjArray(100) ;
616 cpvbranch->SetAddress(&cpvrecpoints);
617 cpvbranch->GetEntry(0);
619 FillRecPointsData(kRPNcpv,cpvrecpoints->GetEntriesFast()) ;
620 cpvrecpoints->Delete();
624 IncEvCountCycleRecPoints();
625 IncEvCountTotalRecPoints();
629 //____________________________________________________________________________
630 void AliPHOSQADataMakerRec::StartOfDetectorCycle()
632 //Detector specific actions at start of cycle