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 **************************************************************************/
17 // Produces the data needed to calculate the quality assurance
18 // All data must be mergeable objects
19 // Handles ESDs and Raws
20 // Histos defined will be used for Raw Data control and monitoring
22 // --- ROOT system ---
23 #include <TClonesArray.h>
30 #include <TParameter.h>
32 // --- Standard library ---
34 // --- AliRoot header files ---
35 #include "AliESDEvent.h"
37 #include "AliCDBManager.h"
38 #include "AliCDBStorage.h"
39 #include "AliCDBEntry.h"
40 #include "AliVZEROQADataMakerRec.h"
41 #include "AliQAChecker.h"
42 #include "AliRawReader.h"
43 #include "AliVZERORawStream.h"
44 #include "AliVZEROdigit.h"
45 #include "AliVZEROReconstructor.h"
49 ClassImp(AliVZEROQADataMakerRec)
51 //____________________________________________________________________________
52 AliVZEROQADataMakerRec::AliVZEROQADataMakerRec() :
53 AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kVZERO), "VZERO Quality Assurance Data Maker"),
60 AliDebug(AliQAv1::GetQADebugLevel(), "Construct VZERO QA Object");
62 for(Int_t i=0; i<64; i++){
66 for(Int_t i=0; i<128; i++){
70 //____________________________________________________________________________
71 AliVZEROQADataMakerRec::AliVZEROQADataMakerRec(const AliVZEROQADataMakerRec& qadm) :
79 SetName((const char*)qadm.GetName()) ;
80 SetTitle((const char*)qadm.GetTitle());
83 //__________________________________________________________________
84 AliVZEROQADataMakerRec& AliVZEROQADataMakerRec::operator = (const AliVZEROQADataMakerRec& qadm )
88 this->~AliVZEROQADataMakerRec();
89 new(this) AliVZEROQADataMakerRec(qadm);
93 //____________________________________________________________________________
94 AliVZEROCalibData* AliVZEROQADataMakerRec::GetCalibData() const
97 AliCDBManager *man = AliCDBManager::Instance();
101 entry = man->Get("VZERO/Calib/Data",fRun);
103 AliWarning("Load of calibration data from default storage failed!");
104 AliWarning("Calibration data will be loaded from local storage ($ALICE_ROOT)");
106 man->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
107 entry = man->Get("VZERO/Calib/Data",fRun);
109 // Retrieval of data in directory VZERO/Calib/Data:
111 AliVZEROCalibData *calibdata = 0;
113 if (entry) calibdata = (AliVZEROCalibData*) entry->GetObject();
114 if (!calibdata) AliFatal("No calibration data from calibration database !");
121 //____________________________________________________________________________
122 void AliVZEROQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
124 // Detector specific actions at end of cycle
125 // Does the QA checking
127 AliQAChecker::Instance()->Run(AliQAv1::kVZERO, task, list) ;
129 for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
130 SetEventSpecie(AliRecoParam::ConvertIndex(specie)) ;
131 if(task == AliQAv1::kRAWS){
132 int nMaxBin = GetRawsData(kPedestalTimeInt0)->GetNbinsY();
133 if(fCurrentCycle%nMaxBin==0) {
134 GetRawsData(kPedestalTimeInt0)->Reset();
135 GetRawsData(kPedestalTimeInt1)->Reset();
136 GetRawsData(kChargeEoITimeInt0)->Reset();
137 GetRawsData(kChargeEoITimeInt1)->Reset();
141 for(Int_t iChannel=0; iChannel<64; iChannel++) {
142 for(Int_t integrator=0;integrator<2;integrator++){
143 sprintf(name,"Ped_%d_%d",iChannel,integrator);
144 hProj = ((TH2I*)GetRawsData((integrator == 0 ? kPedestalCycleInt0 : kPedestalCycleInt1)))->ProjectionY(name,iChannel+1,iChannel+1);
145 ((TH2D*)GetRawsData((integrator == 0 ? kPedestalTimeInt0 : kPedestalTimeInt1)))->Fill((double)iChannel,(double)(fCurrentCycle%nMaxBin),(double)hProj->GetMean());
148 sprintf(name,"Charge_%d_%d",iChannel,integrator);
149 hProj = ((TH2I*)GetRawsData((integrator == 0 ? kChargeEoICycleInt0 : kChargeEoICycleInt1)))->ProjectionY(name,iChannel+1,iChannel+1);
150 ((TH2D*)GetRawsData((integrator == 0 ? kChargeEoITimeInt0 : kChargeEoITimeInt1)))->Fill((double)iChannel,(double)(fCurrentCycle%nMaxBin),hProj->GetMean());
154 } else if (task == AliQAv1::kESDS) {
159 //____________________________________________________________________________
160 void AliVZEROQADataMakerRec::InitESDs()
162 // Creates histograms to control ESDs
164 const Bool_t expert = kTRUE ;
165 const Bool_t image = kTRUE ;
171 h1i = new TH1I("H1I_Cell_Multiplicity_V0A", "Cell Multiplicity in V0A;Multiplicity (Nb of Cell);Counts", 35, 0, 35) ;
172 Add2ESDsList(h1i, kCellMultiV0A, !expert, image) ;
174 h1i = new TH1I("H1I_Cell_Multiplicity_V0C", "Cell Multiplicity in V0;Multiplicity (Nb of Cell);Counts", 35, 0, 35) ;
175 Add2ESDsList(h1i, kCellMultiV0C, !expert, image) ;
177 h1d = new TH1D("H1D_MIP_Multiplicity_V0A", "MIP Multiplicity in V0A;Multiplicity (Nb of MIP);Counts", 1000, 0, 1000) ;
178 Add2ESDsList(h1d, kMIPMultiV0A, !expert, image) ;
180 h1d = new TH1D("H1D_MIP_Multiplicity_V0C", "MIP Multiplicity in V0C;Multiplicity (Nb of MIP);Counts", 1000, 0, 1000) ;
181 Add2ESDsList(h1d, kMIPMultiV0C, !expert, image) ;
183 h2d = new TH2D("H2D_MIP_Multiplicity_Channel", "MIP Multiplicity per Channel;Channel;Multiplicity (Nb of MIP)",64, 0, 64, 100, 0, 100) ;
184 Add2ESDsList(h2d, kMIPMultiChannel, !expert, image) ;
186 h1d = new TH1D("H1D_BBFlag_Counters", "BB Flag Counters;Channel;Counts",64, 0, 64) ;
187 Add2ESDsList(h1d, kBBFlag, !expert, image) ;
189 h1d = new TH1D("H1D_BGFlag_Counters", "BG Flag Counters;Channel;Counts",64, 0, 64) ;
190 Add2ESDsList(h1d, kBGFlag, !expert, image) ;
192 h2d = new TH2D("H2D_Charge_Channel", "ADC Charge per channel;Channel;Charge (ADC counts)",64, 0, 64, 1024, 0, 1024) ;
193 Add2ESDsList(h2d, kChargeChannel, !expert, image) ;
195 h2d = new TH2D("H2D_Time_Channel", "Time per channel;Channel;Time (ns)",64, 0, 64, 820, 0, 410) ;
196 Add2ESDsList(h2d, kTimeChannel, !expert, image) ;
198 h1d = new TH1D("H1D_V0A_Time", "Mean V0A Time;Time (ns);Counts",2048, 0., 409.6);
199 Add2ESDsList(h1d,kESDV0ATime, !expert, image);
201 h1d = new TH1D("H1D_V0C_Time", "Mean V0C Time;Time (ns);Counts",2048, 0., 409.6);
202 Add2ESDsList(h1d,kESDV0CTime, !expert, image);
204 h1d = new TH1D("H1D_Diff_Time", "Diff Time V0A - V0C;Diff Time V0A - V0C (ns);Counts",2*2048, -409.6, 409.6);
205 Add2ESDsList(h1d,kESDDiffTime, !expert, image);
209 //____________________________________________________________________________
210 void AliVZEROQADataMakerRec::InitRaws()
212 // Creates RAW histograms in Raws subdir
214 const Bool_t expert = kTRUE ;
215 const Bool_t saveCorr = kTRUE ;
216 const Bool_t image = kTRUE ;
218 char name[50] , title[100];
219 const Int_t kNintegrator = 2;
221 const Int_t kNTdcTimeBins = 2048;
222 const Int_t kTdcTimeMin = 0;
223 const Int_t kTdcTimeMax = 4096;
224 const Int_t kNTdcWidthBins = 128;
225 const Int_t kTdcWidthMin = 0;
226 const Int_t kTdcWidthMax = 128;
227 const Int_t kNChargeBins = 1024;
228 const Int_t kChargeMin = 0;
229 const Int_t kChargeMax = 1024;
230 const Int_t kNChannelBins = 64;
231 const Int_t kChannelMin = 0;
232 const Int_t kChannelMax = 64;
233 const Int_t kNPedestalBins = 200;
234 const Int_t kPedestalMin = 0;
235 const Int_t kPedestalMax = 200;
236 const Int_t kTimeMin = 0;
237 const Int_t kTimeMax = 100;
238 const Int_t kNMIPBins = 200;
239 const Int_t kMIPMin = 0;
240 const Int_t kMIPMax = 200;
249 // Creation of Cell Multiplicity Histograms
250 h1i = new TH1I("H1I_Multiplicity_V0A", "Cell Multiplicity in V0A;# of Cells;Entries", 35, 0, 35) ;
251 Add2RawsList(h1i,kMultiV0A, !expert, image, saveCorr); iHisto++;
252 h1i = new TH1I("H1I_Multiplicity_V0C", "Cell Multiplicity in V0C;# of Cells;Entries", 35, 0, 35) ;
253 Add2RawsList(h1i,kMultiV0C, !expert, image, saveCorr); iHisto++;
255 // Creation of Total Charge Histograms
256 h1d = new TH1D("H1D_Charge_V0A", "Total Charge in V0A;Charge [ADC counts];Counts", 2048, 0, 32768) ;
257 Add2RawsList(h1d,kChargeV0A, !expert, image, saveCorr); iHisto++;
258 h1d = new TH1D("H1D_Charge_V0C", "Total Charge in V0C;Charge [ADC counts];Counts", 2048, 0, 32768) ;
259 Add2RawsList(h1d,kChargeV0C, !expert, image, saveCorr); iHisto++;
260 h1d = new TH1D("H1D_Charge_V0", "Total Charge in V0;Charge [ADC counts];Counts", 2048, 0, 65536) ;
261 Add2RawsList(h1d,kChargeV0, !expert, image, saveCorr); iHisto++;
263 // Creation of MIP Histograms
264 h1d = new TH1D("H1D_MIP_V0A", "Total MIP in V0A;Charge [MIP];Counts", 2*kNMIPBins,kMIPMin ,32*kMIPMax) ;
265 Add2RawsList(h1d,kRawMIPV0A, !expert, image, saveCorr); iHisto++;
266 h1d = new TH1D("H1D_MIP_V0C", "Total MIP in V0C;Charge [MIP];Counts", 2*kNMIPBins,kMIPMin ,32*kMIPMax) ;
267 Add2RawsList(h1d,kRawMIPV0C, !expert, image, saveCorr); iHisto++;
268 h1d = new TH1D("H1D_MIP_V0", "Total MIP in V0;Charge [MIP];Counts", 2*kNMIPBins,kMIPMin ,32*kMIPMax) ;
269 Add2RawsList(h1d,kRawMIPV0, !expert, image, saveCorr); iHisto++;
270 h2d = new TH2D("H2D_MIP_Channel", "Nb of MIP per channel;Channel;# of Mips", kNChannelBins, kChannelMin, kChannelMax,kNMIPBins,kMIPMin ,kMIPMax) ;
271 Add2RawsList(h2d,kRawMIPChannel, expert, !image, !saveCorr); iHisto++;
274 for(Int_t iInt=0;iInt<kNintegrator;iInt++){
275 // Creation of Pedestal histograms
276 sprintf(name,"H2I_Pedestal_Int%d",iInt);
277 sprintf(title,"Pedestal (Int%d);Pedestal [ADC counts];Counts",iInt);
278 h2i = new TH2I(name, title,kNChannelBins, kChannelMin, kChannelMax,kNPedestalBins,kPedestalMin ,kPedestalMax );
279 Add2RawsList(h2i,(iInt == 0 ? kPedestalInt0 : kPedestalInt1), expert, !image, !saveCorr); iHisto++;
281 // Creation of temporary Pedestal histo used for the mean versus time histogram. This histogram will be reset at the end of each cycle
282 sprintf(name,"H2I_Pedestal_CycleInt%d",iInt);
283 sprintf(title,"One Cycle Pedestal (Int%d);Pedestal [ADC counts];Counts",iInt);
284 h2i = new TH2I(name, title,kNChannelBins, kChannelMin, kChannelMax,kNPedestalBins,kPedestalMin ,kPedestalMax );
285 Add2RawsList(h2i,(iInt == 0 ? kPedestalCycleInt0 : kPedestalCycleInt1), expert, !image, !saveCorr); iHisto++;
287 // Creation of Pedestal versus time graph.
288 sprintf(name,"H2D_Pedestal_Time_Int%d",iInt);
289 sprintf(title,"Pedestal Versus Time (Int%d);Time [ns];Pedestal [ADC counts]",iInt);
290 h2d = new TH2D(name, title,kNChannelBins, kChannelMin, kChannelMax,kTimeMax,kTimeMin ,kTimeMax );
291 Add2RawsList(h2d,(iInt == 0 ? kPedestalTimeInt0 : kPedestalTimeInt1), expert, !image, !saveCorr); iHisto++;
293 // Creation of Charge EoI histograms
294 sprintf(name,"H2I_ChargeEoI_Int%d",iInt);
295 sprintf(title,"Charge EoI (Int%d);Charge [ADC counts];Counts",iInt);
296 h2i = new TH2I(name, title,kNChannelBins, kChannelMin, kChannelMax, kNChargeBins, kChargeMin, kChargeMax);
297 Add2RawsList(h2i,(iInt == 0 ? kChargeEoIInt0 : kChargeEoIInt1), !expert, image, !saveCorr); iHisto++;
299 // Creation of temporary Charge EoI histograms used for the mean versus time histogram. This histogram will be reset at the end of each cycle
300 sprintf(name,"H2I_ChargeEoI_CycleInt%d",iInt);
301 sprintf(title,"One Cycle Charge EoI (Int%d);Charge [ADC counts];Counts",iInt);
302 h2i = new TH2I(name, title,kNChannelBins, kChannelMin, kChannelMax, kNChargeBins, kChargeMin, kChargeMax);
303 Add2RawsList(h2i,(iInt == 0 ? kChargeEoICycleInt0 : kChargeEoICycleInt1), expert, !image, !saveCorr); iHisto++;
305 // Creation of Charge EoI versus time graphs
306 sprintf(name,"H2D_ChargeEoI_Time_Int%d",iInt);
307 sprintf(title,"Charge EoI Versus Time (Int%d);Time [ns];Charge [ADC counts]",iInt);
308 h2d = new TH2D(name, title,kNChannelBins, kChannelMin, kChannelMax,kTimeMax,kTimeMin ,kTimeMax );
309 Add2RawsList(h2d,(iInt == 0 ? kChargeEoITimeInt0 : kChargeEoITimeInt1), expert, !image, !saveCorr); iHisto++;
311 sprintf(name,"H2I_ChargeEoI_BB_Int%d",iInt);
312 sprintf(title,"Charge EoI w/ BB Flag (Int%d);??;Charge [ADC counts]",iInt);
313 h2i = new TH2I(name, title,kNChannelBins, kChannelMin, kChannelMax, kNChargeBins, kChargeMin, kChargeMax);
314 Add2RawsList(h2i,(iInt == 0 ? kChargeEoIBBInt0 : kChargeEoIBBInt1), expert, !image, !saveCorr); iHisto++;
316 sprintf(name,"H2I_ChargeEoI_BG_Int%d",iInt);
317 sprintf(title,"Charge EoI w/ BG Flag (Int%d);??;Charge [ADC counts]",iInt);
318 h2i = new TH2I(name, title,kNChannelBins, kChannelMin, kChannelMax, kNChargeBins, kChargeMin, kChargeMax);
319 Add2RawsList(h2i,(iInt == 0 ? kChargeEoIBGInt0: kChargeEoIBGInt1), expert, !image, !saveCorr); iHisto++;
321 // Creation of Charge versus LHC Clock histograms
322 sprintf(name,"H2D_ChargeVsClock_Int%d",iInt);
323 sprintf(title,"Charge Versus LHC-Clock (Int%d);Tine [ns];Charge [ADC counts]",iInt);
324 h2d = new TH2D(name, title,kNChannelBins, kChannelMin, kChannelMax,21, -10.5, 10.5 );
325 Add2RawsList(h2d,(iInt == 0 ? kChargeVsClockInt0 : kChargeVsClockInt1 ), expert, !image, !saveCorr); iHisto++;
327 // Creation of Minimum Bias Charge histograms
328 for(Int_t iBB=0;iBB<2;iBB++){
329 for(Int_t iBG=0;iBG<2;iBG++){
330 sprintf(name,"H2I_ChargeMB_BB%d_BG%d_Int%d",iBB,iBG,iInt);
331 sprintf(title,"MB Charge (BB=%d, BG=%d, Int=%d);Charge [ADC counts];Counts",iBB,iBG,iInt);
332 h2i = new TH2I(name, title,kNChannelBins, kChannelMin, kChannelMax,kNChargeBins, kChargeMin, kChargeMax);
336 if(iBG==0) idx = kChargeMBBB0BG0Int0;
337 else idx = kChargeMBBB0BG1Int0;
339 if(iBG==0) idx = kChargeMBBB1BG0Int0;
340 else idx = kChargeMBBB1BG1Int0;
344 if(iBG==0) idx = kChargeMBBB0BG0Int1;
345 else idx = kChargeMBBB0BG1Int1;
347 if(iBG==0) idx = kChargeMBBB1BG0Int1;
348 else idx = kChargeMBBB1BG1Int1;
351 Add2RawsList(h2i,idx, expert, !image, !saveCorr); iHisto++;
357 // Creation of Time histograms
358 sprintf(name,"H2I_Width");
359 sprintf(title,"HPTDC Width;Width [ns];Counts");
360 h2i = new TH2I(name, title,kNChannelBins, kChannelMin, kChannelMax, kNTdcWidthBins, kTdcWidthMin, kTdcWidthMax);
361 Add2RawsList(h2i,kWidth, expert, !image, !saveCorr); iHisto++;
363 sprintf(name,"H2I_Width_BB");
364 sprintf(title,"HPTDC Width w/ BB Flag condition;??;Width [ns]");
365 h2i = new TH2I(name, title,kNChannelBins, kChannelMin, kChannelMax, kNTdcWidthBins, kTdcWidthMin, kTdcWidthMax);
366 Add2RawsList(h2i,kWidthBB, expert, !image, !saveCorr); iHisto++;
368 sprintf(name,"H2I_Width_BG");
369 sprintf(title,"HPTDC Width w/ BG Flag condition??;Width [ns]");
370 h2i = new TH2I(name, title,kNChannelBins, kChannelMin, kChannelMax, kNTdcWidthBins, kTdcWidthMin, kTdcWidthMax);
371 Add2RawsList(h2i,kWidthBG, expert, !image, !saveCorr); iHisto++;
373 sprintf(name,"H2I_HPTDCTime");
374 sprintf(title,"HPTDC Time;??;Width [ns]");
375 h2i = new TH2I(name, title,kNChannelBins, kChannelMin, kChannelMax, kNTdcTimeBins, kTdcTimeMin, kTdcTimeMax);
376 Add2RawsList(h2i,kHPTDCTime, !expert, image, !saveCorr); iHisto++;
378 sprintf(name,"H2I_HPTDCTime_BB");
379 sprintf(title,"HPTDC Time w/ BB Flag condition;??;Width [ns]");
380 h2i = new TH2I(name, title,kNChannelBins, kChannelMin, kChannelMax, kNTdcTimeBins, kTdcTimeMin, kTdcTimeMax);
381 Add2RawsList(h2i,kHPTDCTimeBB, expert, !image, !saveCorr); iHisto++;
383 sprintf(name,"H2I_HPTDCTime_BG");
384 sprintf(title,"HPTDC Time w/ BG Flag condition;??;Width [ns]");
385 h2i = new TH2I(name, title,kNChannelBins, kChannelMin, kChannelMax, kNTdcTimeBins, kTdcTimeMin, kTdcTimeMax);
386 Add2RawsList(h2i,kHPTDCTimeBG, expert, !image, !saveCorr); iHisto++;
388 sprintf(name,"H1D_V0A_Time");
389 sprintf(title,"V0A Time;Time [ns];Counts");
390 h1d = new TH1D(name, title,kNTdcTimeBins, kTdcTimeMin/10, kTdcTimeMax/10);
391 Add2RawsList(h1d,kV0ATime, !expert, !image, saveCorr); iHisto++;
393 sprintf(name,"H1D_V0C_Time");
394 sprintf(title,"V0C Time;Time [ns];Counts");
395 h1d = new TH1D(name, title,kNTdcTimeBins, kTdcTimeMin/10, kTdcTimeMax/10);
396 Add2RawsList(h1d,kV0CTime, !expert, !image, saveCorr); iHisto++;
398 sprintf(name,"H1D_Diff_Time");
399 sprintf(title,"Diff V0A-V0C Time;Time [ns];Counts");
400 h1d = new TH1D(name, title,2*kNTdcTimeBins, -kTdcTimeMax/10, kTdcTimeMax/10);
401 Add2RawsList(h1d,kDiffTime, !expert, !image, saveCorr); iHisto++;
403 // Creation of Flag versus LHC Clock histograms
404 sprintf(name,"H2D_BBFlagVsClock");
405 sprintf(title,"BB-Flags Versus LHC-Clock;Time [ns];??");
406 h2d = new TH2D(name, title,kNChannelBins, kChannelMin, kChannelMax,21, -10.5, 10.5 );
407 Add2RawsList(h2d,kBBFlagVsClock, expert, !image, !saveCorr); iHisto++;
409 sprintf(name,"H2D_BGFlagVsClock");
410 sprintf(title,"BG-Flags Versus LHC-Clock;Time [ns];??");
411 h2d = new TH2D(name, title,kNChannelBins, kChannelMin, kChannelMax,21, -10.5, 10.5 );
412 Add2RawsList(h2d,kBGFlagVsClock, expert, !image, !saveCorr); iHisto++;
414 AliDebug(AliQAv1::GetQADebugLevel(), Form("%d Histograms has been added to the Raws List",iHisto));
417 //____________________________________________________________________________
418 void AliVZEROQADataMakerRec::InitDigits()
420 // create Digits histograms in Digits subdir
421 const Bool_t expert = kTRUE ;
422 const Bool_t image = kTRUE ;
430 // create Digits histograms in Digits subdir
431 TH1I * h0 = new TH1I("hDigitMultiplicity", "Digits multiplicity distribution in VZERO;# of Digits;Entries", 100, 0, 99) ;
433 Add2DigitsList(h0, 0, !expert, image) ;
435 for (Int_t i=0; i<64; i++)
437 sprintf(TDCname, "hDigitTDC%d", i);
438 sprintf(texte,"Digit TDC in cell %d; TDC value;Entries",i);
439 fhDigTDC[i] = new TH1I(TDCname,texte,300,0.,149.);
441 sprintf(ADCname,"hDigitADC%d",i);
442 sprintf(texte,"Digit ADC in cell %d;ADC value;Entries",i);
443 fhDigADC[i]= new TH1I(ADCname,texte,1024,0.,1023.);
445 Add2DigitsList(fhDigTDC[i],i+1, !expert, image);
446 Add2DigitsList(fhDigADC[i],i+1+64, !expert, image);
450 //____________________________________________________________________________
451 void AliVZEROQADataMakerRec::MakeDigits(TClonesArray * digits)
453 // makes data from Digits
455 GetDigitsData(0)->Fill(digits->GetEntriesFast()) ;
457 AliVZEROdigit *VZERODigit ;
458 while ( (VZERODigit = dynamic_cast<AliVZEROdigit *>(next())) ) {
459 Int_t PMNumber = VZERODigit->PMNumber();
460 GetDigitsData(PMNumber +1)->Fill( VZERODigit->Time()) ; // in 100 of picoseconds
461 GetDigitsData(PMNumber +1+64)->Fill( VZERODigit->ADC()) ;
466 //____________________________________________________________________________
467 void AliVZEROQADataMakerRec::MakeDigits(TTree *digitTree)
469 // makes data from Digit Tree
471 TClonesArray * digits = new TClonesArray("AliVZEROdigit", 1000) ;
473 TBranch * branch = digitTree->GetBranch("VZERODigit") ;
475 AliWarning("VZERO branch in Digit Tree not found") ;
477 branch->SetAddress(&digits) ;
478 branch->GetEntry(0) ;
484 //____________________________________________________________________________
485 void AliVZEROQADataMakerRec::MakeESDs(AliESDEvent * esd)
487 // Creates QA data from ESDs
489 UInt_t eventType = esd->GetEventType();
493 AliESDVZERO *esdVZERO=esd->GetVZEROData();
495 if (!esdVZERO) break;
497 GetESDsData(kCellMultiV0A)->Fill(esdVZERO->GetNbPMV0A());
498 GetESDsData(kCellMultiV0C)->Fill(esdVZERO->GetNbPMV0C());
499 GetESDsData(kMIPMultiV0A)->Fill(esdVZERO->GetMTotV0A());
500 GetESDsData(kMIPMultiV0C)->Fill(esdVZERO->GetMTotV0C());
502 Float_t timeV0A = 0., timeV0C = 0., diffTime;
503 Int_t iTimeV0A = 0, iTimeV0C = 0;
505 for(Int_t i=0;i<64;i++) {
506 GetESDsData(kMIPMultiChannel)->Fill((Float_t) i,(Float_t) esdVZERO->GetMultiplicity(i));
507 GetESDsData(kChargeChannel)->Fill((Float_t) i,(Float_t) esdVZERO->GetAdc(i));
508 if(esdVZERO->GetBBFlag(i)) GetESDsData(kBBFlag)->Fill((Float_t) i);
509 if(esdVZERO->GetBGFlag(i)) GetESDsData(kBGFlag)->Fill((Float_t) i);
511 Float_t time = (Float_t) esdVZERO->GetTime(i)/10.; //Convert in ns: 1 TDC channel = 100ps
512 GetESDsData(kTimeChannel)->Fill((Float_t) i,time);
524 if(iTimeV0A>0) timeV0A /= iTimeV0A;
526 if(iTimeV0C>0) timeV0C /= iTimeV0C;
528 if(timeV0A<0. || timeV0C<0.) diffTime = -10000.;
529 else diffTime = timeV0A - timeV0C;
531 GetESDsData(kESDV0ATime)->Fill(timeV0A);
532 GetESDsData(kESDV0CTime)->Fill(timeV0C);
533 GetESDsData(kESDDiffTime)->Fill(diffTime);
540 //____________________________________________________________________________
541 void AliVZEROQADataMakerRec::MakeRaws(AliRawReader* rawReader)
543 // Fills histograms with Raws, computes average ADC values dynamically (pedestal subtracted)
546 AliVZERORawStream* rawStream = new AliVZERORawStream(rawReader);
549 eventTypeType eventType = rawReader->GetType();
553 Double_t timeV0A =0., timeV0C = 0.;
554 UInt_t itimeV0A=0, itimeV0C=0;
555 Double_t chargeV0A=0., chargeV0C=0.;
556 Double_t mipV0A=0., mipV0C=0.;
558 Double_t diffTime=-100000.;
569 Int_t MBCharge, charge;
573 for(Int_t iChannel=0; iChannel<64; iChannel++) { // BEGIN : Loop over channels
575 offlineCh = rawStream->GetOfflineChannel(iChannel);
577 // Fill Pedestal histograms
579 for(Int_t j=15; j<21; j++) {
580 if((rawStream->GetBGFlag(iChannel,j) || rawStream->GetBBFlag(iChannel,j))) iFlag++;
583 if(iFlag == 0){ //No Flag found
584 for(Int_t j=15; j<21; j++){
585 pedestal=rawStream->GetPedestal(iChannel, j);
586 integrator = rawStream->GetIntegratorFlag(iChannel, j);
588 GetRawsData((integrator == 0 ? kPedestalInt0 : kPedestalInt1))->Fill(offlineCh,pedestal);
589 GetRawsData((integrator == 0 ? kPedestalCycleInt0 : kPedestalCycleInt1))->Fill(offlineCh,pedestal);
593 // Fill Charge EoI histograms
595 // Look for the maximum in the LHC clock train
599 for(Int_t iEvent=0; iEvent<21; iEvent++){
600 iCharge = rawStream->GetPedestal(iChannel,iEvent);
605 } // End of maximum searching procedure
607 integrator = rawStream->GetIntegratorFlag(iChannel,iClock);
608 BBFlag = rawStream->GetBBFlag(iChannel, iClock);
609 BGFlag = rawStream->GetBGFlag(iChannel,iClock );
611 GetRawsData((integrator == 0 ? kChargeEoIInt0 : kChargeEoIInt1))->Fill(offlineCh,charge);
612 if(BBFlag) GetRawsData((integrator == 0 ? kChargeEoIBBInt0 : kChargeEoIBBInt1))->Fill(offlineCh,charge);
613 if(BGFlag) GetRawsData((integrator == 0 ? kChargeEoIBGInt0 : kChargeEoIBGInt1))->Fill(offlineCh,charge);
615 hProj = ((TH2I*)GetRawsData((integrator == 0 ? kPedestalInt0 : kPedestalInt1)))->ProjectionY("",offlineCh+1,offlineCh+1);
616 Double_t ped = hProj->GetMean();
617 Double_t sigma = hProj->GetRMS();
620 Double_t chargeEoI = charge - ped;
622 // Calculation of the number of MIP
623 Double_t mipEoI = chargeEoI * fCalibData->GetMIPperADC(offlineCh);
626 if(charge<1023 && chargeEoI > 5.*sigma){
627 ((TH2I*)GetRawsData((integrator == 0 ? kChargeEoICycleInt0 : kChargeEoICycleInt1)))->Fill(offlineCh,chargeEoI);
628 ((TH2D*)GetRawsData(kRawMIPChannel))->Fill(offlineCh,mipEoI);
631 chargeV0C += chargeEoI;
635 chargeV0A += chargeEoI;
640 // Fill Charge Minimum Bias Histograms
643 for(Int_t iBunch=0; iBunch<10; iBunch++){
644 integrator = rawStream->GetIntMBFlag(iChannel, iBunch);
645 BBFlag = rawStream->GetBBMBFlag(iChannel, iBunch);
646 BGFlag = rawStream->GetBGMBFlag(iChannel, iBunch);
647 MBCharge = rawStream->GetChargeMB(iChannel, iBunch);
651 if(BGFlag==0) idx = kChargeMBBB0BG0Int0;
652 else idx = kChargeMBBB0BG1Int0;
654 if(BGFlag==0) idx = kChargeMBBB1BG0Int0;
655 else idx = kChargeMBBB1BG1Int0;
659 if(BGFlag==0) idx = kChargeMBBB0BG0Int1;
660 else idx = kChargeMBBB0BG1Int1;
662 if(BGFlag==0) idx = kChargeMBBB1BG0Int1;
663 else idx = kChargeMBBB1BG1Int1;
666 GetRawsData(idx)->Fill(offlineCh,MBCharge);
669 // Fill HPTDC Time Histograms
671 BBFlag = rawStream->GetBBFlag(iChannel, 10);
672 BGFlag = rawStream->GetBGFlag(iChannel, 10);
673 time = rawStream->GetTime(iChannel);
674 width = rawStream->GetWidth(iChannel);
685 GetRawsData(kHPTDCTime)->Fill(offlineCh,time);
686 GetRawsData(kWidth)->Fill(offlineCh,width);
688 GetRawsData(kHPTDCTimeBB)->Fill(offlineCh,time);
689 GetRawsData(kWidthBB)->Fill(offlineCh,width);
692 GetRawsData(kHPTDCTimeBG)->Fill(offlineCh,time);
693 GetRawsData(kWidthBG)->Fill(offlineCh,width);
696 // Fill Flag and Charge Versus LHC-Clock histograms
698 for(Int_t iEvent=0; iEvent<21; iEvent++){
699 charge = rawStream->GetPedestal(iChannel,iEvent);
700 integrator = rawStream->GetIntegratorFlag(iChannel,iEvent);
701 BBFlag = rawStream->GetBBFlag(iChannel, iEvent);
702 BGFlag = rawStream->GetBGFlag(iChannel,iEvent );
704 ((TH2*) GetRawsData((integrator == 0 ? kChargeVsClockInt0 : kChargeVsClockInt1 )))->Fill(offlineCh,(float)iEvent-10,(float)charge);
705 ((TH2*) GetRawsData(kBBFlagVsClock))->Fill(offlineCh,(float)iEvent-10,(float)BBFlag);
706 ((TH2*) GetRawsData(kBGFlagVsClock))->Fill(offlineCh,(float)iEvent-10,(float)BGFlag);
709 }// END of Loop over channels
711 if(itimeV0A>0) timeV0A /= (itimeV0A * 10); // itimeV0A Channels and divide by 10 to have the result in ns because 1 TDC Channel = 100 ps
713 if(itimeV0C>0) timeV0C /= (itimeV0C * 10);
715 if(timeV0A<0. || timeV0C<0.) diffTime = -10000.;
716 else diffTime = timeV0A - timeV0C;
718 GetRawsData(kV0ATime)->Fill(timeV0A);
719 GetRawsData(kV0CTime)->Fill(timeV0C);
720 GetRawsData(kDiffTime)->Fill(diffTime);
722 GetRawsData(kMultiV0A)->Fill(mulV0A);
723 GetRawsData(kMultiV0C)->Fill(mulV0C);
725 GetRawsData(kChargeV0A)->Fill(chargeV0A);
726 GetRawsData(kChargeV0C)->Fill(chargeV0C);
727 GetRawsData(kChargeV0)->Fill(chargeV0A + chargeV0C);
729 GetRawsData(kRawMIPV0A)->Fill(mipV0A);
730 GetRawsData(kRawMIPV0C)->Fill(mipV0C);
731 GetRawsData(kRawMIPV0)->Fill(mipV0A + mipV0C);
734 } // END of SWITCH : EVENT TYPE
737 TParameter<double> * p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kMultiV0A)->GetName()))) ;
738 if (p) p->SetVal((double)mulV0A) ;
740 p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kMultiV0C)->GetName()))) ;
741 if (p) p->SetVal((double)mulV0C) ;
743 p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kChargeV0A)->GetName()))) ;
744 if (p) p->SetVal((double)chargeV0A) ;
746 p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kChargeV0C)->GetName()))) ;
747 if (p) p->SetVal((double)chargeV0C) ;
749 p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kChargeV0)->GetName()))) ;
750 if (p) p->SetVal((double)(chargeV0A + chargeV0C)) ;
752 p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kRawMIPV0A)->GetName()))) ;
753 if (p) p->SetVal((double)mipV0A) ;
755 p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kRawMIPV0C)->GetName()))) ;
756 if (p) p->SetVal((double)mipV0C) ;
758 p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kRawMIPV0)->GetName()))) ;
759 if (p) p->SetVal((double)(mipV0A + mipV0C)) ;
761 p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kV0ATime)->GetName()))) ;
762 if (p) p->SetVal((double)timeV0A) ;
764 p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kV0CTime)->GetName()))) ;
765 if (p) p->SetVal((double)timeV0C) ;
767 p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kDiffTime)->GetName()))) ;
768 if (p) p->SetVal((double)diffTime) ;
770 delete rawStream; rawStream = 0x0;
775 //____________________________________________________________________________
776 void AliVZEROQADataMakerRec::StartOfDetectorCycle()
778 // Detector specific actions at start of cycle
780 // Reset of the histogram used - to have the trend versus time -
782 fCalibData = GetCalibData();
785 h = GetRawsData(kPedestalCycleInt0);
787 h = GetRawsData(kPedestalCycleInt1);
789 h = GetRawsData(kChargeEoICycleInt0);
791 h = GetRawsData(kChargeEoICycleInt1);