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 **************************************************************************/
18 Produces the data needed to calculate the quality assurance.
19 All data must be mergeable objects.
22 // --- ROOT system ---
23 #include <TClonesArray.h>
28 #include <TParameter.h>
30 // --- Standard library ---
32 // --- AliRoot header files ---
33 #include "AliESDEvent.h"
35 #include "AliCDBManager.h"
36 #include "AliCDBStorage.h"
37 #include "AliCDBEntry.h"
38 #include "AliVZEROQADataMakerRec.h"
39 #include "AliQAChecker.h"
40 #include "AliRawReader.h"
41 #include "AliVZERORawStream.h"
42 #include "AliVZEROReconstructor.h"
45 ClassImp(AliVZEROQADataMakerRec)
47 //____________________________________________________________________________
48 AliVZEROQADataMakerRec::AliVZEROQADataMakerRec() :
49 AliQADataMakerRec(AliQA::GetDetName(AliQA::kVZERO), "VZERO Quality Assurance Data Maker"),
50 fCalibData(GetCalibData()),
56 for(Int_t i=0; i<64; i++){
60 for(Int_t i=0; i<128; i++){
64 //____________________________________________________________________________
65 AliVZEROQADataMakerRec::AliVZEROQADataMakerRec(const AliVZEROQADataMakerRec& qadm) :
67 fCalibData(GetCalibData()),
72 SetName((const char*)qadm.GetName()) ;
73 SetTitle((const char*)qadm.GetTitle());
76 //__________________________________________________________________
77 AliVZEROQADataMakerRec& AliVZEROQADataMakerRec::operator = (const AliVZEROQADataMakerRec& qadm )
80 this->~AliVZEROQADataMakerRec();
81 new(this) AliVZEROQADataMakerRec(qadm);
85 //____________________________________________________________________________
86 AliVZEROCalibData* AliVZEROQADataMakerRec::GetCalibData() const
89 AliCDBManager *man = AliCDBManager::Instance();
93 entry = man->Get("VZERO/Calib/Data");
95 // Retrieval of data in directory VZERO/Calib/Data:
97 AliVZEROCalibData *calibdata = 0;
99 if (entry) calibdata = (AliVZEROCalibData*) entry->GetObject();
100 if (!calibdata) AliFatal("No calibration data from calibration database !");
105 //____________________________________________________________________________
106 void AliVZEROQADataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX_t task, TObjArray * list)
108 //Detector specific actions at end of cycle
109 // do the QA checking
110 AliQAChecker::Instance()->Run(AliQA::kVZERO, task, list) ;
113 //____________________________________________________________________________
114 void AliVZEROQADataMakerRec::InitESDs()
116 //Create histograms to control ESD
118 TH1I * h1 = new TH1I("hVZERONbPMA", "Number of PMs fired in V0A", 80, 0, 80);
120 Add2ESDsList(h1, 0) ;
122 TH1I * h2 = new TH1I("hVZERONbPMC", "Number of PMs fired in V0C", 80, 0, 80);
124 Add2ESDsList(h2, 1) ;
126 TH1I * h3 = new TH1I("hVZEROMultA", "Multiplicity in V0A", 50, 0, 50) ;
128 Add2ESDsList(h3, 2) ;
130 TH1I * h4 = new TH1I("hVZEROMultC", "Multiplicity in V0C", 50, 0, 50) ;
132 Add2ESDsList(h4, 3) ;
134 TH2F * h5 = new TH2F("fVzeroMult", "Vzero multiplicity",
135 64, -0.5, 63.5,1000, -0.5, 99.5);
136 h5->GetXaxis()->SetTitle("Vzero PMT");
137 h5->GetYaxis()->SetTitle("Multiplicity");
139 Add2ESDsList(h5, 4) ;
140 TH1F * h6 = new TH1F("fBBA","BB Vzero A", 32, -0.5,31.5);
142 Add2ESDsList(h6, 5) ;
143 TH1F * h7 = new TH1F("fBGA","BG Vzero A", 32, -0.5,31.5);
145 Add2ESDsList(h7, 6) ;
146 TH1F * h8 = new TH1F("fBBC","BB Vzero C", 32, -0.5,31.5);
148 Add2ESDsList(h8, 7) ;
149 TH1F * h9 = new TH1F("fBGC","BG Vzero C", 32, -0.5,31.5);
151 Add2ESDsList(h9, 8) ;
153 TH2F *h10 = new TH2F("fVzeroAdc", "Vzero Adc",
154 64, -0.5, 63.5,1024, -0.5, 1023.5);
155 h10->GetXaxis()->SetTitle("Vzero PMT");
156 h10->GetYaxis()->SetTitle("Adc counts");
158 Add2ESDsList(h10, 9);
160 TH2F *h11 = new TH2F("fVzeroTime", "Vzero Time",
161 64, -0.5, 63.5,300, -0.5, 149.5);
162 h11->GetXaxis()->SetTitle("Vzero PMT");
163 h11->GetYaxis()->SetTitle("Time [100 ps]");
165 Add2ESDsList(h11,10);
169 //____________________________________________________________________________
170 void AliVZEROQADataMakerRec::InitRaws()
172 // create Raws histograms in Raws subdir
179 Bool_t expert = kTRUE ;
180 Bool_t saveCorr = kTRUE ;
182 TH2I *h0 = new TH2I("hCellADCMap0","ADC vs Cell for EVEN Integrator", 70, 0, 70, 512, 0, 1024);
184 Add2RawsList(h0,0, !expert, !saveCorr) ;
185 TH2I *h1 = new TH2I("hCellADCMap1","ADC vs Cell for ODD Integrator", 70, 0, 70, 512, 0, 1024);
187 Add2RawsList(h1,1, !expert, !saveCorr) ;
189 TH2F *hMeanADC0 = new TH2F("hCellMeanADCMap0","Mean ADC vs cell for EVEN integrator",70,-0.5,69.5,512,-0.5,511.5);
190 Add2RawsList(hMeanADC0,2, !expert, !saveCorr);
191 TH2F *hMeanADC1 = new TH2F("hCellMeanADCMap1","Mean ADC vs cell for ODD integrator",70,-0.5,69.5,512,-0.5,511.5);
192 Add2RawsList(hMeanADC1,3, !expert, !saveCorr);
194 TH1I *hMulV0A = new TH1I("hMulV0A", "Multiplicity in V0A", 40, 0, 40) ;
195 Add2RawsList(hMulV0A,4, !expert, saveCorr);
196 TH1I *hMulV0C = new TH1I("hMulV0C", "Multiplicity in V0C", 40, 0, 40) ;
197 Add2RawsList(hMulV0C,5, !expert, saveCorr);
199 for (Int_t i=0; i<64; i++)
201 sprintf(ADCname,"hRaw0ADC%d",i);
202 sprintf(texte,"Raw ADC in cell %d for even integrator",i);
203 hRawADC0[i]= new TH1I(ADCname,texte,1024,0,1024);
204 Add2RawsList(hRawADC0[i],i+6, expert, !saveCorr);
206 sprintf(ADCname,"hRaw1ADC%d",i);
207 sprintf(texte,"Raw ADC in cell %d for odd integrator",i);
208 hRawADC1[i]= new TH1I(ADCname,texte,1024,0,1024);
209 Add2RawsList(hRawADC1[i],i+6+64, expert, !saveCorr);
213 //____________________________________________________________________________
214 void AliVZEROQADataMakerRec::MakeESDs(AliESDEvent * esd)
216 // make QA data from ESDs
218 AliESDVZERO *esdVZERO=esd->GetVZEROData();
221 GetESDsData(0)->Fill(esdVZERO->GetNbPMV0A());
222 GetESDsData(1)->Fill(esdVZERO->GetNbPMV0C());
223 GetESDsData(2)->Fill(esdVZERO->GetMTotV0A());
224 GetESDsData(3)->Fill(esdVZERO->GetMTotV0C());
225 for(Int_t i=0;i<64;i++) {
226 GetESDsData(4)->Fill((Float_t) i,(Float_t) esdVZERO->GetMultiplicity(i));
227 GetESDsData(9)->Fill((Float_t) i,(Float_t) esdVZERO->GetAdc(i));
228 GetESDsData(10)->Fill((Float_t) i,(Float_t) esdVZERO->GetTime(i));
230 for(Int_t i=0;i<32;i++) {
231 if(esdVZERO->BBTriggerV0A(i))
232 GetESDsData(5)->Fill((Float_t) i);
233 if(esdVZERO->BGTriggerV0A(i))
234 GetESDsData(6)->Fill((Float_t) i);
235 if(esdVZERO->BBTriggerV0C(i))
236 GetESDsData(7)->Fill((Float_t) i);
237 if(esdVZERO->BGTriggerV0C(i))
238 GetESDsData(8)->Fill((Float_t) i);
244 //____________________________________________________________________________
245 void AliVZEROQADataMakerRec::MakeRaws(AliRawReader* rawReader)
247 //Fill histograms with Raws, computes average ADC values dynamically (pedestal subtracted)
250 AliVZERORawStream* rawStream = new AliVZERORawStream(rawReader);
253 Float_t ChargeEoI, Threshold; // for pedestal subtraction
255 // for(Int_t i=0; i<128; i++) { printf(" i = %d pedestal = %f sigma = %f \n\n",
256 // i, fCalibData->GetPedestal(i),fCalibData->GetSigma(i) );}
258 GetRawsData(2)->Reset(); // to keep only the last value of the ADC average
259 GetRawsData(3)->Reset();
264 Int_t thresholV0A = 50 ;
265 Int_t thresholV0C = 50 ;
267 for(Int_t j=0; j<64; j++) {
268 Int_t i = rawStream->GetOfflineChannel(j); // Convert Online to Offline channel number
269 if(!rawStream->GetIntegratorFlag(i,10))
271 // even integrator - fills index6 to 69
272 GetRawsData(0)->Fill(i,rawStream->GetADC(i)) ;
273 GetRawsData(i+6)->Fill(rawStream->GetADC(i)) ;
274 ChargeEoI = (float)rawStream->GetADC(i) - fCalibData->GetPedestal(i);
275 Threshold = 3.0 * fCalibData->GetSigma(i);
276 if(rawStream->GetADC(i)<1023 && ChargeEoI > Threshold) {
277 fADC_Mean[i] = ((fADC_Mean[i]*fEven[i]) + rawStream->GetADC(i)) / (fEven[i]+1);
278 GetRawsData(2)->Fill((Float_t) i, fADC_Mean[i]);
281 else if(rawStream->GetIntegratorFlag(i,10))
283 // odd integrator - fills index 70 to 133
284 GetRawsData(1)->Fill(i,rawStream->GetADC(i)) ;
285 GetRawsData(i+6+64)->Fill(rawStream->GetADC(i)) ;
286 ChargeEoI = (float)rawStream->GetADC(i) - fCalibData->GetPedestal(i+64) ;
287 Threshold = 3.0 * fCalibData->GetSigma(i+64);
288 if(rawStream->GetADC(i)<1023 && ChargeEoI > Threshold) {
289 fADC_Mean[i+64] = ((fADC_Mean[i+64]*fOdd[i]) + rawStream->GetADC(i)) / (fOdd[i]+1);
290 GetRawsData(3)->Fill((Float_t) i, fADC_Mean[i+64]);
294 if (rawStream->GetADC(i) > thresholV0C)
297 if (rawStream->GetADC(i) > thresholV0A)
303 GetRawsData(4)->Fill(mulV0A) ;
304 TParameter<double> * p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQA::GetTaskName(AliQA::kRAWS).Data(), GetRawsData(4)->GetName()))) ;
306 GetRawsData(5)->Fill(mulV0C) ;
307 p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQA::GetTaskName(AliQA::kRAWS).Data(), GetRawsData(5)->GetName()))) ;
311 //____________________________________________________________________________
312 void AliVZEROQADataMakerRec::StartOfDetectorCycle()
314 // Detector specific actions at start of cycle