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 // Alla.Maevskaya@cern.ch
24 // --- ROOT system ---
25 #include <TClonesArray.h>
29 #include <TDirectory.h>
31 // --- Standard library ---
33 // --- AliRoot header files ---
34 #include "AliESDEvent.h"
36 #include "AliT0digit.h"
38 #include "AliT0RecPoint.h"
39 #include "AliT0QADataMakerRec.h"
40 #include "AliQAChecker.h"
41 #include "AliT0RawReader.h"
42 #include "AliT0RecoParam.h"
43 #include "THnSparse.h"
45 #include "Riostream.h"
46 ClassImp(AliT0QADataMakerRec)
48 //____________________________________________________________________________
49 AliT0QADataMakerRec::AliT0QADataMakerRec() :
50 AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kT0),
51 "T0 Quality Assurance Data Maker"),
56 for (Int_t i=0; i<6; i++) {
63 for (Int_t i=0; i<24; i++)
74 // for(Int_t ic=0; ic<24; ic++)
75 // fMeans[ic] = new TH1F(Form("CFD1minCFD%d",ic+1),"CFD-CFD",100,-250,250);
79 //____________________________________________________________________________
80 AliT0QADataMakerRec::AliT0QADataMakerRec(const AliT0QADataMakerRec& qadm) :
86 SetName((const char*)qadm.GetName()) ;
87 SetTitle((const char*)qadm.GetTitle());
90 //__________________________________________________________________
91 AliT0QADataMakerRec& AliT0QADataMakerRec::operator = (const AliT0QADataMakerRec& qadm )
94 this->~AliT0QADataMakerRec();
95 new(this) AliT0QADataMakerRec(qadm);
98 //__________________________________________________________________
99 AliT0QADataMakerRec::~AliT0QADataMakerRec()
103 //____________________________________________________________________________
104 void AliT0QADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
106 //Detector specific actions at end of cycle
107 // do the QA checking
108 AliQAChecker::Instance()->Run(AliQAv1::kT0, task, list) ;
110 for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
111 if (! IsValidEventSpecie(specie, list))
113 SetEventSpecie(AliRecoParam::ConvertIndex(specie)) ;
114 if ( task == AliQAv1::kRAWS ) {
115 GetRawsData(0)->SetLabelSize(0.02);
116 const Char_t *triggers[6] = {"mean", "vertex","ORA","ORC","central","semi-central"};
117 for (Int_t itr=0; itr<6; itr++) {
119 fTrEffCal[itr] = Float_t (fNumTriggersCal[itr])/Float_t (fnEventCal);
121 fTrEffPhys[itr] = Float_t (fNumTriggers[itr])/Float_t (fnEventPhys);
123 GetRawsData(169+250)->Fill(triggers[itr], fTrEffCal[itr]);
124 GetRawsData(169+250)->SetBinContent(itr+1, fTrEffCal[itr]);
125 GetRawsData(169)->Fill(triggers[itr], fTrEffPhys[itr]);
126 GetRawsData(169)->SetBinContent(itr+1, fTrEffPhys[itr]);
129 for(Int_t ik=0; ik<24; ik++)
132 if ( fnEventCal>0) effic = Float_t(feffC[ik])/Float_t(fnEventCal);
133 GetRawsData(207+250)->SetBinContent(ik+1,effic) ;
135 if ( fnEventCal>0) effic = Float_t(feffA[ik])/Float_t(fnEventCal);
136 GetRawsData(208+250)->SetBinContent(ik+1,effic );
138 if ( fnEventCal>0) effic = Float_t(feffqtc[ik])/Float_t(fnEventCal);
139 GetRawsData(209+250)->SetBinContent(ik+1, effic);
142 if ( fnEventPhys>0) effic = Float_t(feffPhysC[ik])/Float_t(fnEventPhys);
143 GetRawsData(207)->SetBinContent(ik+1, effic);
153 //____________________________________________________________________________
154 void AliT0QADataMakerRec::StartOfDetectorCycle()
156 //Detector specific actions at start of cycle
160 //____________________________________________________________________________
161 void AliT0QADataMakerRec::InitRaws()
163 // create Raw histograms in Raw subdir
164 const Bool_t expert = kTRUE ;
165 const Bool_t saveCorr = kTRUE ;
166 const Bool_t image = kTRUE ;
171 for (Int_t i=0; i<500; i++){
177 TString timename, ampname, qtcname, ledname;
178 TString timeCalname, ampCalname, ledCalname, qtcCalname;
179 TString qt1name, qt0name, qt1Calname, qt0Calname;
182 TH1F* fhRefPoint = new TH1F("hRefPoint","Ref Point", 10000, 0 ,50000);
183 Add2RawsList( fhRefPoint,0, expert, !image, !saveCorr);
185 TH1F* fhRefPointcal = new TH1F("hRefPointcal","Ref Point laser", 5000, 0 ,20000);
186 Add2RawsList( fhRefPointcal,250, expert, !image, !saveCorr);
189 TH1F * fhRawLEDamp[24];
190 TH1F *fhRawQTC[24]; TH1F * fhRawLED[24];
191 TH1F *fhRawCFDcal[24]; TH1F * fhRawLEDampcal[24];
192 TH1F *fhRawQTCcal[24]; TH1F * fhRawLEDcal[24];
193 TH1F *fhRawQT0cal[24]; TH1F *fhRawQT1cal[24];
194 TH1F *fhRawQT1[24]; TH1F *fhRawQT0[24];
195 TH1F* fhRawNhits[24];
197 for (Int_t i=0; i<24; i++)
202 qt0name = "hRawQT0_";
203 qt1name = "hRawQT1_";
204 ampname = "hRawLEDminCFD";
213 fhRawCFD[i] = new TH1F(timename.Data(), Form("%s;CFD [#channels];Counts", timename.Data()),Int_t((high[i+1]-low[i+1])/4),low[i+1],high[i+1]);
214 Add2RawsList( fhRawCFD[i],i+1, expert, !image, !saveCorr);
215 fhRawLED[i] = new TH1F(ledname.Data(), Form("%s;LED[#channels];Counts", ledname.Data()),Int_t((high[i+25]-low[i+25])/4),low[i+25],high[i+25]);
216 Add2RawsList( fhRawLED[i],i+25, expert, !image, !saveCorr);
217 fhRawLEDamp[i] = new TH1F(ampname.Data(), Form("%s;LED-CFD [#channels];Counts", ampname.Data()),1000,0,1000);
218 Add2RawsList( fhRawLEDamp[i],i+49, expert, !image, !saveCorr);
219 fhRawQTC[i] = new TH1F(qtcname.Data(), Form("%s;QTC[#channels];Counts", qtcname.Data()),10000,0,10000);
220 Add2RawsList( fhRawQTC[i],i+73, expert, !image, !saveCorr);
221 fhRawQT1[i] = new TH1F(qt1name.Data(), Form("%s;QT1[#channels];Counts", qt1name.Data()),Int_t((high[97+i]-low[97+i])/4),low[97+i],high[97+i]);
222 Add2RawsList( fhRawQT1[i],97+i, expert, !image, !saveCorr);
223 fhRawQT0[i] = new TH1F(qt0name.Data(), Form("%s;QT0[#channels];Counts", qt0name.Data()),Int_t((high[121+i]-low[121+i])/4),low[121+i],high[121+i]);
224 Add2RawsList( fhRawQT0[i],121+i, expert, !image, !saveCorr);
226 fhRawNhits[i] = new TH1F(nhits.Data(), Form("%s;#Hits;Events", nhits.Data()),20, 0, 20);
227 Add2RawsList( fhRawNhits[i],176+i, expert, !image, !saveCorr);
230 timeCalname ="hRawCFDcal";
231 ledCalname = "hRawLEDcal";
232 ampCalname = "hRawLEDminCFDcal";
233 qtcCalname = "hRawQTCcal";
234 qt0Calname = "hRawQT0cal";
235 qt1Calname = "hRawQT1cal";
242 fhRawCFDcal[i] = new TH1F(timeCalname.Data(), Form("LASER: %s;Time ;Counts", timeCalname.Data()),Int_t((high[251+i]-low[251+i])/4),low[251+i],high[251+i]);
243 Add2RawsList( fhRawCFDcal[i],251+i, expert, !image, !saveCorr);
245 fhRawLEDcal[i] = new TH1F(ledCalname.Data(), Form("LASER: %s;Time ;Counts", ledCalname.Data()),Int_t((high[251+i+24]-low[251+i+24])/4),low[251+i+24],high[251+i+24]);
246 Add2RawsList( fhRawLEDcal[i],251+i+24, expert, !image, !saveCorr);
248 fhRawLEDampcal[i] = new TH1F(ampCalname.Data(), Form("LASER: %s;Amplitude [ADC counts];Counts", ampCalname.Data()),1000,0,1000);
249 Add2RawsList( fhRawLEDampcal[i],251+i+48, expert, !image, !saveCorr);
251 fhRawQTCcal[i] = new TH1F(qtcCalname.Data(), Form("LASER: %s;QTC[#channels]; ;Counts",qtcCalname.Data()),10000,0,10000);
252 Add2RawsList( fhRawQTCcal[i],251+i+72, expert, !image, !saveCorr);
254 fhRawQT0cal[i] = new TH1F(qt0Calname.Data(), Form("LASER: %s;QT0[#channels] ;Counts",qt0Calname.Data()),Int_t((high[251+96+i]-low[251+96+i])/4),low[251+96+i],high[251+96+i]);
255 Add2RawsList( fhRawQT0cal[i],251+96+i, expert, !image, !saveCorr);
257 fhRawQT1cal[i] = new TH1F(qt1Calname.Data(), Form("LASER: %s;QT1[#channels] ;Counts",qt1Calname.Data()),Int_t((high[i+251+120]-low[i+251+120])/4),low[i+251+120],high[i+251+120]);
258 Add2RawsList( fhRawQT1cal[i],i+251+120, expert, !image, !saveCorr);
263 TH1F* fhRawTrigger = new TH1F("hRawTrigger"," phys triggers;Trigger ;Counts",6,0,6);
264 Add2RawsList(fhRawTrigger ,169, expert, image, !saveCorr);
266 TH1F* fhRawMean = new TH1F("hRawMean","online mean signal, physics event;",Int_t((high[170]-low[170])/4),low[170],high[170]);
267 Add2RawsList( fhRawMean,170, expert, !image, !saveCorr);
269 TH1F* fhRawVertex = new TH1F("hRawVertex","online vertex signal; counts",Int_t((high[171]-low[171])/4),low[171],high[171]);
270 Add2RawsList( fhRawVertex,171, expert, image, !saveCorr);
272 TH1F* fhRawORA = new TH1F("hRawORA","online OR A; counts",Int_t((high[172]-low[172])/4),low[172],high[172]);
273 Add2RawsList( fhRawORA,172, expert, !image, !saveCorr);
274 TH1F* fhRawORC = new TH1F("hRawORC","online OR C;counts",Int_t(( high[173]-low[173])/4),low[173],high[173]);
275 Add2RawsList( fhRawORC,173, expert, !image, !saveCorr);
276 TH1F* fhMultCentr = new TH1F("hMultCentr","online trigger Central;counts ",Int_t(( high[174]-low[174])/4),low[174],high[174]);
277 Add2RawsList( fhMultCentr,174, expert, !image, !saveCorr);
278 TH1F* fhMultSeCentr = new TH1F("hMultSemiCentr","online trigger SemiCentral;counts ",Int_t(( high[175]-low[175])/4),low[175],high[175]);
279 Add2RawsList( fhMultSeCentr,175, expert, !image, !saveCorr);
282 TH1F* fhRawTriggerCal = new TH1F("hRawTriggerCal"," laser triggers",6,0,6);
283 Add2RawsList(fhRawTriggerCal ,169+250 , !expert, image, !saveCorr);
285 TH1F* fhRawMeanCal = new TH1F("hRawMeanCal","online mean signal, calibration event",Int_t((high[170+250]-low[170+250])/4),low[170+250],high[170+250]);
286 Add2RawsList( fhRawMeanCal,170+250, expert, !image, !saveCorr);
288 TH1F* fhRawVertexCal = new TH1F("hRawVertexCal","online vertex signal, calibration event ",Int_t((high[171+250]-low[171+250])/4),low[171+250],high[171+250] );
289 Add2RawsList( fhRawVertexCal,171+250, expert, !image, !saveCorr);
292 TH1F* fhRawORAcal = new TH1F("hRawORAcal","laser OR A; counts",Int_t((high[172+250]-low[172+250])/4),low[172+250],high[172+250]);
293 Add2RawsList( fhRawORAcal,172+250, expert, !image, !saveCorr );
296 TH1F* fhRawORCcal = new TH1F("hRawORCcal","laserOR C;counts ",Int_t(( high[173]-low[173])/4),low[173],high[173]);
297 Add2RawsList( fhRawORCcal,173+250, expert, !image, !saveCorr);
299 TH1F* fhMultCentrcal = new TH1F("hMultCentrcal","laser trigger Central;counts ",Int_t(( high[174]-low[174])/4),low[174],high[174]);
300 Add2RawsList( fhMultCentrcal,174+250, expert, !image, !saveCorr);
301 TH1F* fhMultSeCentrcal = new TH1F("hMultSemiCentrcal","laser trigger SemiCentral;counts ",Int_t(( high[175]-low[175])/4),low[175],high[175]);
302 Add2RawsList( fhMultSeCentrcal,175+250, expert, !image, !saveCorr);
304 //multiplicity trigger
306 TH1F* fhMultAcal = new TH1F("hMultAcal","laser: full mulltiplicity;Multiplicity A side;Entries",Int_t((high[201]-low[201])/4),low[201],high[201]);
307 Add2RawsList( fhMultAcal,201+250, expert, !image, !saveCorr );
308 TH1F* fhMultAScal = new TH1F("hMultASemical","laser:full multiplicity with semi-central trigger A side;Multiplicity;Entries",
309 Int_t((high[202]-low[202])/4),low[202],high[202] );
310 Add2RawsList( fhMultAScal,202+250, expert, !image, !saveCorr);
311 TH1F* fhMultACcal = new TH1F("hMultACentrcal","laser:full multiplicity with central trigger A side;Multiplicity;Entries",
312 Int_t((high[203]-low[203])/4),low[203],high[203]);
313 Add2RawsList( fhMultACcal,203+250, expert, !image, !saveCorr);
315 TH1F* fhMultA = new TH1F("hMultA","full mulltiplicity A side;Multiplicity;Entries", Int_t((high[201]-low[201])/4) ,low[201],high[201]);
316 Add2RawsList( fhMultA,201, expert, image, !saveCorr );
318 TH1F* fhMultAS = new TH1F("hMultASemi","full multiplicity with semi-central trigger A side ;Multiplicity;Entries",
319 Int_t((high[202]-low[202])/4),low[202],high[202] );
320 Add2RawsList( fhMultAS, 202, expert, !image, !saveCorr);
321 TH1F* fhMultAC = new TH1F("hMultACentr","full multiplicity with central trigger;Multiplicity;Entries",
322 Int_t((high[203]-low[203])/4),low[203],high[203]);
323 Add2RawsList( fhMultAC, 203, expert, !image, !saveCorr);
327 TH1F* fhMultCcal = new TH1F("hMultCcal","laser:full mulltiplicity C side;Multiplicity;Entries",Int_t((high[204]-low[204])/4),low[204],high[204]);
328 Add2RawsList( fhMultCcal,204+250, expert, !image, !saveCorr );
329 TH1F* fhMultCScal = new TH1F("hMultCSemical","laser:full multiplicity with semi-central trigger C side;Multiplicity;Entries",
330 Int_t((high[205]-low[205])/4),low[205],high[205] );
331 Add2RawsList( fhMultCScal,205+250, expert, !image, !saveCorr);
332 TH1F* fhMultCCcal = new TH1F("hMultCCentrcal","laser:full multiplicity with central trigger C side;Multiplicity;Entries",
333 Int_t((high[206]-low[206])/4),low[206],high[206]);
334 Add2RawsList( fhMultCCcal,206+250, expert, !image, !saveCorr);
336 TH1F* fhMultC = new TH1F("hMultC","full mulltiplicity C side;Multiplicity;Entries", Int_t(high[204]-low[204]/4) ,low[204],high[204]);
337 Add2RawsList( fhMultC,204, expert, image, !saveCorr );
338 TH1F* fhMultCS = new TH1F("hMultCSemi","full multiplicity with semi-central trigger C side;Multiplicity;Entries",
339 Int_t((high[205]-low[205])/4),low[205],high[205] );
340 Add2RawsList( fhMultCS,205, expert, !image, !saveCorr);
341 TH1F* fhMultCC = new TH1F("hMultCentr","full multiplicity with central trigger C side;Multiplicity;Entries",
342 Int_t((high[206]-low[206])/4),low[206],high[206]);
343 Add2RawsList( fhMultCC,206, expert, !image, !saveCorr);
347 TH1F* fhEffCFD = new TH1F("hEffCFDcal","CFD efficiecy laser ;#PMT; #CFD counts/nEvents",24, 0 ,24);
348 Add2RawsList( fhEffCFD,207+250, !expert, image, !saveCorr);
350 TH1F* fhCFDeffpsys= new TH1F("fhCFDeffpsys"," CFD efficiency; #PMT; #CFD counts/nEvents",24, 0 ,24);
351 // fhCFDeffpsys->SetMaximum(2);
352 Add2RawsList( fhCFDeffpsys, 207, expert, image, !saveCorr);
354 TH1F* fhEffLED = new TH1F("hEffLEDcal","LEDefficiecy; #PMT; #LED counts/nEvent",24, 0 ,24);
355 Add2RawsList( fhEffLED,208+250, !expert, image, !saveCorr);
357 TH1F* fhEffQTC = new TH1F("hEffQTCcal","QTC efficiecy; #PMT; QTC efficiency%s;",24, 0 ,24);
358 Add2RawsList( fhEffQTC,209+250, !expert, image, !saveCorr);
361 TH2F* fhCFD = new TH2F("hCFD","CFD phys; #PMT; CFD {#channnels}",25, 0 ,25,Int_t((high[210]-low[210])/4),low[210],high[210]);
362 fhCFD->SetOption("COLZ");
363 Add2RawsList( fhCFD,210, expert, image, !saveCorr);
365 TH2F* fhLED = new TH2F("hLED","LED phys; #PMT; LED [#channnels]",25, 0 ,25,Int_t((high[211]-low[211])/4),low[211],high[211]);
366 fhLED->SetOption("COLZ");
367 Add2RawsList( fhLED,211, expert, image, !saveCorr);
369 TH2F* fhQTC = new TH2F("hQTC","QTC phys; #PMT; QTC [#channnels]",25, 0 ,25,Int_t( high[212]-low[212]),low[212],high[212]);
370 fhQTC->SetOption("COLZ");
371 Add2RawsList( fhQTC,212, expert, image, !saveCorr);
373 TH2F* fhCFDcal = new TH2F("hCFDcal","CFD laser; #PMT; CFD {#channnels}",25, 0 ,25,Int_t((high[210]-low[210])/4),low[210],high[210]);
374 fhCFDcal->SetOption("COLZ");
375 Add2RawsList( fhCFDcal,210+250, expert, image, !saveCorr);
378 TH2F* fhLEDcal = new TH2F("hLEDcal","LED laser; #PMT; LED [#channnels]",25, 0 ,25,Int_t((high[211]-low[211])/4),low[211],high[211]);
379 fhLEDcal->SetOption("COLZ");
380 Add2RawsList( fhLEDcal,211+250, expert, image, !saveCorr);
382 TH2F* fhQTCcal = new TH2F("hQTCcal","QTC laser; #PMT; QTC [#channnels]",25, 0 ,25,Int_t( high[212]-low[212]),low[212],high[212]);
383 fhQTCcal->SetOption("COLZ");
384 Add2RawsList( fhQTCcal,212+250, expert, image, !saveCorr);
387 TH1F* fhNumPMTA= new TH1F("hNumPMTA","number of PMT hitted per event",13, 0 ,13);
388 Add2RawsList(fhNumPMTA ,213, expert, image, !saveCorr);
390 TH1F* fhNumPMTC= new TH1F("hNumPMTC","number of PMT hitted per event",13, 0 ,13);
391 Add2RawsList(fhNumPMTC ,214, expert, image, !saveCorr);
393 TH1F* fhHitsOrA= new TH1F("fhHitsOrA","T0_OR A hit multiplicitie",20, 0 ,20);
394 Add2RawsList( fhHitsOrA,215, expert, !image, !saveCorr);
396 TH1F* fhHitsOrC= new TH1F("fhHitsOrC","T0_OR C hit multiplicitie",20, 0 ,20);
397 Add2RawsList(fhHitsOrC ,216, expert, !image, !saveCorr);
400 TH1F* fhOrCminOrA= new TH1F("fhOrCminOrA","T0_OR C - T0_OR A",10000,-5000,5000);
401 Add2RawsList( fhOrCminOrA,219, expert, !image, !saveCorr);
403 TH1F* fhOrCminOrAcal= new TH1F("fhOrCminOrAcal","T0_OR C - T0_OR A",10000,-5000,5000);
404 Add2RawsList( fhOrCminOrAcal,219+250, expert, !image, !saveCorr);
406 TH1F* fhOrCminOrATvdcOn= new TH1F("fhOrCminOrATvdcOn","T0_OR C - T0_OR A TVDC on",10000,-5000,5000);
407 Add2RawsList( fhOrCminOrATvdcOn,217, expert, !image, !saveCorr);
409 TH1F* fhOrCminOrATvdcOncal= new TH1F("fhOrCminOrATvdcOncal","T0_OR C - T0_OR A TVDC on laser",10000,-5000,5000);
410 Add2RawsList( fhOrCminOrATvdcOncal,217+250, expert, !image, !saveCorr);
412 TH1F* fhOrCminOrATvdcOff= new TH1F("fhOrCminOrATvdcOff","T0_OR C - T0_OR A TVDC off",10000,-5000,5000);
413 Add2RawsList( fhOrCminOrATvdcOff,218, expert, !image, !saveCorr);
416 TH1F* fhOrCminOrATvdcOffcal= new TH1F("fhOrCminOrATvdcOffcal","T0_OR C - T0_OR ATVDC off laser",10000,-5000,5000);
417 Add2RawsList( fhOrCminOrATvdcOffcal,218+250, expert, !image, !saveCorr);
419 //satellite & beam background
420 TH2F* fhBeam = new TH2F("fhBeam", " Mean vs Vertex ", 120, -30, 30, 120, -30, 30);
421 Add2RawsList( fhBeam,220, !expert, image, !saveCorr);
422 TH2F* fhBeamTVDCon = new TH2F("fhBeamTVDCon", " Mean vs Vertex TVDC on ",50, -5, 5, 50, -5, 5) ;
423 Add2RawsList( fhBeamTVDCon,221, expert, image, !saveCorr);
424 TH2F* fhBeamTVDCoff = new TH2F("fhBeamTVDCoff", " Mean vs Vertex TVDC off", 120, -30, 30, 120, -30, 30);
425 Add2RawsList( fhBeamTVDCoff,222, expert, image, !saveCorr);
426 TH1F* fhMean = new TH1F("fhMean", " (T0A+T0C)/2 ", 200, -2000, 2000);
427 Add2RawsList( fhMean,223, !expert, image, !saveCorr);
430 const Char_t *triggers[6] = {"mean", "vertex","ORA","ORC","central","semi-central"};
431 for (Int_t itr=0; itr<6; itr++) {
432 GetRawsData(169)->Fill(triggers[itr], fNumTriggersCal[itr]);
433 GetRawsData(169)->SetBinContent(itr+1, fNumTriggersCal[itr]);
434 GetRawsData(169+250)->Fill(triggers[itr], fNumTriggers[itr]);
435 GetRawsData(169+250)->SetBinContent(itr+1, fNumTriggers[itr]);
441 //____________________________________________________________________________
442 void AliT0QADataMakerRec::InitDigits()
444 // create Digits histograms in Digits subdir
445 const Bool_t expert = kTRUE ;
446 const Bool_t image = kTRUE ;
448 TH2F * fhDigCFD = new TH2F("fhDigCFD", " CFD digits; #PMT; CFD digits[#channels]",25,-0.5,24.5,100,0,1000);
449 fhDigCFD->SetOption("COLZ");
450 Add2DigitsList( fhDigCFD,0, !expert, image);
451 TH2F *fhDigLEDamp = new TH2F("fhDigLEDamp", " LED-CFD digits; #PMT; LED-CFD amplitude ",25,-0.5,24.5,100,100,1000);
452 fhDigLEDamp->SetOption("COLZ");
453 Add2DigitsList( fhDigLEDamp,1, !expert, !image);
454 TH2F * fhDigQTC = new TH2F("fhDigQTC", " QTC digits; #PMT; QTC amplitude",25,-0.5,24.5,100,100,10000);
455 fhDigQTC->SetOption("COLZ");
456 Add2DigitsList( fhDigQTC,2, !expert, !image);
461 //____________________________________________________________________________
463 void AliT0QADataMakerRec::InitRecPoints()
465 // create cluster histograms in RecPoint subdir
466 const Bool_t expert = kTRUE ;
467 const Bool_t image = kTRUE ;
469 TH2F* fhRecCFD = new TH2F("hRecCFD"," CFD time;#PMT; CFD Time [ns];",24, 0 ,24,
471 fhRecCFD->SetOption("COLZ");
472 Add2RecPointsList ( fhRecCFD,0, !expert, image);
474 TH2F* fhRecAmpDiff = new TH2F("hRecAmpDiff"," LED-CFD min QTC amplitude;#PMT; difference [MIPs];",
475 24, 0 ,24, 200,-10,10);
476 fhRecAmpDiff->SetOption("COLZ");
477 Add2RecPointsList (fhRecAmpDiff, 1, !expert, image);
479 TH1F *fhMean = new TH1F("hMean","online - rec mean;online - rec mean[#channels];",2000, -1000, 1000);
480 Add2RecPointsList ( fhMean,2, !expert, image);
484 //____________________________________________________________________________
485 void AliT0QADataMakerRec::InitESDs()
487 //create ESDs histograms in ESDs subdir
488 const Bool_t expert = kTRUE ;
489 const Bool_t image = kTRUE ;
491 TH1F *fhESDMean = new TH1F("hESDmean"," ESD mean; mean time[%channels]",1000, -5, 5);
492 Add2ESDsList(fhESDMean, 0, expert, !image) ;
493 TH1F * fhESDVertex = new TH1F("hESDvertex","ESDvertex; vertex[cm];",82,-30,30);
494 Add2ESDsList(fhESDVertex, 1, expert, !image) ;
496 TH1F * fhESDResolution = new TH1F("hESDResolution","(T0A-T0C)/2 corrected by SPD vertex; ns",800,-2,2);
497 Add2ESDsList(fhESDResolution, 2, !expert, image) ;
501 //____________________________________________________________________________
502 void AliT0QADataMakerRec::MakeRaws( AliRawReader* rawReader)
505 for(Int_t i=0; i<24; i++) time[i] = 0;
508 //fills QA histos for RAW
510 // Int_t refPointParam = GetRecoParam()->GetRefPoint();
512 Int_t refPointParam = 0;
514 AliT0RawReader *start = new AliT0RawReader(rawReader);
517 AliDebug(AliQAv1::GetQADebugLevel(),Form(" no raw data found!!"));
520 UInt_t type =rawReader->GetType();
521 if (type == 8){ shift=250; fnEventCal++;}
522 if (type == 7){ shift=0; fnEventPhys++;}
523 // if (type == 7){ shift=1; fnEventPhys++;}
524 Int_t allData[110][5];
525 for (Int_t i0=0; i0<110; i0++)
527 for (Int_t j0=0; j0<5; j0++) allData[i0][j0]=0;
529 for (Int_t i=0; i<107; i++)
530 for (Int_t iHit=0; iHit<5; iHit++)
531 allData[i][iHit]= start->GetData(i,iHit);
533 if ( allData[0][0] > 0 && (type == 7))
534 GetRawsData(0) -> Fill( allData[0][0]);
535 if ( allData[0][0] > 0 && (type == 8))
536 GetRawsData(250) -> Fill( allData[0][0]);
538 refpoint = allData[refPointParam][0];
539 if (refPointParam < 0 ) refpoint=0;
540 if (refPointParam == 0 ) refpoint = allData[0][0] - 5000;
542 Int_t sideshift, sideshiftqtc;
545 for (Int_t ik = 0; ik<24; ik++)
550 if(allData[ik+sideshift][0]>0 && type == 7 ) numPmtC++;
553 if(allData[ik+45][0]>0 && type == 7 ) numPmtA++;
560 for (Int_t iHt=0; iHt<5; iHt++){
562 if(allData[ik+sideshift][iHt]>0) {
563 GetRawsData(shift+ik+1) -> Fill(allData[ik+sideshift][iHt]);
564 GetRawsData(210+shift)->Fill(ik+1, allData[ik+sideshift][iHt]);
565 if(type == 8 ) feffC[ik]++;
566 AliDebug(50,Form("%i CFD %i data %s",ik, ik+sideshift, GetRawsData(shift+ik+1)->GetName()));
574 if(allData[ik+12+sideshift][iHt] > 0) {
575 GetRawsData(shift+ik+24+1)-> Fill(allData[ik+12+sideshift][iHt]);
576 GetRawsData(211+shift)->Fill(ik+1, allData[ik+12+sideshift][iHt]);
577 AliDebug(50,Form("%i LED %i data %s",ik, ik+12+sideshift, GetRawsData(shift+ik+1+24)->GetName()));
584 if(allData[ik+12+sideshift][iHt] > 0 && allData[ik+sideshift][iHt] >0 )
585 GetRawsData(shift+ik+48+1)->
586 Fill(allData[ik+12+sideshift][iHt]-allData[ik+sideshift][iHt]);
589 if(allData[2*ik+sideshiftqtc+24][iHt] > 0 &&
590 allData[2*ik+sideshiftqtc+25][iHt] > 0) {
591 GetRawsData(shift+ik+72+1)->
592 Fill(allData[2*ik+sideshiftqtc+24][iHt]-allData[2*ik+sideshiftqtc+25][iHt]);
593 GetRawsData(212+shift)->Fill(ik+1, allData[2*ik+sideshiftqtc+24][iHt]-allData[2*ik+sideshiftqtc+25][iHt]);
594 if(type == 8) feffqtc[ik]++;
595 AliDebug(50,Form("%i QTC %i data %s",ik, 2*ik+sideshiftqtc+24, GetRawsData(shift+ik+1+72)->GetName()));
598 if(allData[2*ik+sideshiftqtc+24][iHt] > 0) {
599 AliDebug(50,Form("%i QT0 %i data %s",ik, 2*ik+sideshiftqtc+24, GetRawsData(shift+ik+1+96)->GetName()));
600 GetRawsData(shift+ik+96+1)->Fill(allData[2*ik+sideshiftqtc+24][iHt]);
602 if(allData[2*ik+sideshiftqtc+25][iHt] > 0) {
603 AliDebug(50,Form("%i QT0 %i data %s",ik, 2*ik+sideshiftqtc+25, GetRawsData(shift+ik+1+120)->GetName()));
604 GetRawsData(shift+ik+120+1)->Fill(allData[2*ik+sideshiftqtc+25][iHt]);
609 GetRawsData(ik+176)->Fill(nhitsPMT);
610 GetRawsData(213)->Fill(numPmtC);
611 GetRawsData(214)->Fill(numPmtA);
615 Int_t trChannel[6] = {49,50,51,52,55,56};
616 Float_t ch2cm = 24.4*0.029979;
619 for (Int_t iHt=0; iHt<5; iHt++) {
621 //orA-orC phys tvdc 1
622 if((allData[51][iHt]>0 && allData[52][iHt]>0) && allData[50][iHt]>0)
624 AliDebug(10,Form("orA-orC phys tvdc 1 %i data %s", 217+shift, GetRawsData(shift+217)->GetName()));
626 GetRawsData(217+shift)->Fill((allData[52][iHt]-allData[51][iHt])*ch2cm);
627 // GetRawsData(345) ->Fill((allData[51][iHt]+allData[52][iHt])/2.);
629 //orA-orC phys tvdc 0
630 if((allData[51][iHt]>0 && allData[52][iHt]>0) && allData[50][iHt]<=0)
632 AliDebug(10,Form("orA-orC phys tvdc 0 %i data %s", 218+shift, GetRawsData(shift+218)->GetName()));
634 GetRawsData(218+shift)->Fill((allData[52][iHt]-allData[51][iHt])*ch2cm);
636 if(allData[51][iHt]>0 && allData[52][iHt]>0) {
637 AliDebug(50,Form("orA-orC phys tvdc all %i data %s", 219+shift, GetRawsData(shift+219)->GetName()));
638 GetRawsData(219+shift)->Fill((allData[52][iHt]-allData[51][iHt])*ch2cm);
640 for (Int_t itr=0; itr<6; itr++) {
641 if (allData[trChannel[itr]][iHt] >0) {
642 if(type == 7 )fNumTriggers[itr]++;
643 if(type == 8 )fNumTriggersCal[itr]++;
644 AliDebug(50,Form(" triggers %i data %s", 170+itr+shift, GetRawsData(170+itr+shift)->GetName()));
646 GetRawsData(170+itr+shift)->Fill(allData[trChannel[itr]][iHt]);
650 if(type == 7) if(allData[51][iHt] >0) nhitsOrA++;
651 if(type == 7)if(allData[52][iHt] >0) nhitsOrC++;
653 //mult trigger signals phys
655 if(allData[53][iHt]>0 && allData[54][iHt]>0) {
656 AliDebug(50,Form(" mpdA %i data %s", 201+shift, GetRawsData(201+shift)->GetName()));
658 GetRawsData(201+shift)->Fill(allData[53][iHt]-allData[54][iHt]);
659 if(allData[56][iHt]>0) GetRawsData(202+shift)->Fill(allData[53][iHt]-allData[54][iHt]);
660 if(allData[55][iHt]>0) GetRawsData(203+shift)->Fill(allData[53][iHt]-allData[54][iHt]);
664 if(allData[105][iHt]>0 && allData[106][iHt]>0) {
665 AliDebug(50,Form(" mpdC %i data %s", 204+shift, GetRawsData(204+shift)->GetName()));
667 GetRawsData(204+shift)->Fill(allData[105][iHt]-allData[106][iHt]);
668 if(allData[56][iHt]>0) GetRawsData(205+shift)->Fill(allData[105][iHt]-allData[106][iHt]);
669 if(allData[55][iHt]>0) GetRawsData(206+shift)->Fill(allData[105][iHt]-allData[106][iHt]);
673 GetRawsData(215)->Fill(nhitsOrA);
674 GetRawsData(216)->Fill(nhitsOrC);
678 if(fnEventPhys == 1000) {
679 for (Int_t ik=0; ik<24; ik++)
681 GetRawsData(ik+1)->GetXaxis()->SetRangeUser(2000, 3000);
682 int maxBin = GetRawsData(ik+1) ->GetMaximumBin();
683 double meanEstimate = GetRawsData(ik+1)->GetBinCenter( maxBin);
684 TF1* fit= new TF1("fit","gaus", meanEstimate - 40, meanEstimate + 40);
685 fit->SetParameters (GetRawsData(ik+1)->GetBinContent(maxBin), meanEstimate, 80);
686 GetRawsData(ik+1) -> Fit("fit","RQ","Q", meanEstimate-40, meanEstimate+40);
687 fMeans[ik]= (Int_t) fit->GetParameter(1);
688 GetRawsData(ik+1)->GetXaxis()->SetRangeUser(0, 30000);
693 if( fnEventPhys > 1000) {
694 Int_t besttimeA=9999999;
695 Int_t besttimeC=9999999;
696 for (Int_t ipmt=0; ipmt<12; ipmt++){
697 if(allData[ipmt+1][0] > 1 ) {
698 // time[ipmt] = allData[ipmt+1][0] - Int_t(GetRawsData(1)->GetMean());
699 time[ipmt] = allData[ipmt+1][0] - fMeans[ipmt];
700 if(TMath::Abs(time[ipmt]) < TMath::Abs(besttimeC))
701 besttimeC=time[ipmt]; //timeC
704 for ( Int_t ipmt=12; ipmt<24; ipmt++){
705 if(allData[ipmt+45][0] > 0) {
706 time[ipmt] = allData[ipmt+45][0] - fMeans[ipmt] ;
707 if(TMath::Abs(time[ipmt]) < TMath::Abs(besttimeA) );
708 besttimeA=time[ipmt]; //timeA
711 if(besttimeA<99999 &&besttimeC< 99999) {
712 Float_t t0 = 24.4 * (Float_t( besttimeA+besttimeC)/2. );
713 Float_t ver = 24.4 * Float_t( besttimeA-besttimeC)/2.;
714 GetRawsData(220)->Fill(0.001*ver, 0.001*(t0));
715 if(allData[50][0] > 0) GetRawsData(221)->Fill(0.001*ver, 0.001*(t0));
716 if(allData[50][0] <= 0) GetRawsData(222)->Fill(0.001*ver, 0.001*(t0));
717 GetRawsData(223) ->Fill(t0);
731 //____________________________________________________________________________
732 void AliT0QADataMakerRec::MakeDigits( TTree *digitsTree)
734 //fills QA histos for Digits
736 TArrayI *digCFD = new TArrayI(24);
737 TArrayI *digLED = new TArrayI(24);
738 TArrayI *digQT0 = new TArrayI(24);
739 TArrayI *digQT1 = new TArrayI(24);
742 TBranch *brDigits=digitsTree->GetBranch("T0");
743 AliT0digit *fDigits = new AliT0digit() ;
745 brDigits->SetAddress(&fDigits);
747 AliError(Form("EXEC Branch T0 digits not found"));
750 digitsTree->GetEvent(0);
751 digitsTree->GetEntry(0);
752 brDigits->GetEntry(0);
753 fDigits->GetTimeCFD(*digCFD);
754 fDigits->GetTimeLED(*digLED);
755 fDigits->GetQT0(*digQT0);
756 fDigits->GetQT1(*digQT1);
757 refpoint = fDigits->RefPoint();
758 for (Int_t i=0; i<24; i++)
760 if (digCFD->At(i)>0) {
761 Int_t cfd=digCFD->At(i)- refpoint;
762 GetDigitsData(0) ->Fill(i,cfd);
763 GetDigitsData(1) -> Fill(i, (digLED->At(i) - digCFD->At(i)));
764 GetDigitsData(2) -> Fill(i, (digQT1->At(i) - digQT0->At(i)));
775 //____________________________________________________________________________
776 void AliT0QADataMakerRec::MakeRecPoints(TTree * clustersTree)
778 //fills QA histos for clusters
780 AliT0RecPoint* frecpoints= new AliT0RecPoint ();
782 AliError(":MakeRecPoints >> no recpoints found");
785 TBranch *brRec =clustersTree ->GetBranch("T0");
787 brRec->SetAddress(&frecpoints);
789 AliError(Form("EXEC Branch T0 rec not found "));
795 for ( Int_t i=0; i<24; i++) {
797 GetRecPointsData(0) -> Fill(i, frecpoints -> GetTime(i) - frecpoints -> GetTime(0));
799 GetRecPointsData(0) -> Fill(i, frecpoints -> GetTime(i) - frecpoints -> GetTime(12));
800 GetRecPointsData(1) -> Fill( i, frecpoints -> GetAmp(i) - frecpoints->AmpLED(i));
802 Double_t mmm=frecpoints->GetOnlineMean()- frecpoints->GetMeanTime();
803 GetRecPointsData(2) ->Fill(mmm);
807 //____________________________________________________________________________
808 void AliT0QADataMakerRec::MakeESDs(AliESDEvent * esd)
810 //fills QA histos for ESD
812 const Double32_t *mean;
813 mean = esd->GetT0TOF();
814 Double32_t t0time= 0.001*mean[0];
815 Double32_t orA= 0.001*mean[1];
816 Double32_t orC=0.001* mean[2];
818 if (t0time<99) GetESDsData(0) -> Fill(t0time);
819 if( esd->GetT0zVertex() <99) GetESDsData(1)-> Fill(esd->GetT0zVertex());
820 if( orA<99 && orC<99) GetESDsData(2)-> Fill((orA-orC)/2.);
823 //____________________________________________________________________________
826 void AliT0QADataMakerRec::GetMeanAndSigma(TH1F* hist, Float_t &mean, Float_t &sigma)
829 const double window = 3.; //fit window
831 double meanEstimate, sigmaEstimate;
833 maxBin = hist->GetMaximumBin(); //position of maximum
834 meanEstimate = hist->GetBinCenter( maxBin); // mean of gaussian sitting in maximum
835 sigmaEstimate = hist->GetRMS();
836 TF1* fit= new TF1("fit","gaus", meanEstimate - window*sigmaEstimate, meanEstimate + window*sigmaEstimate);
837 fit->SetParameters(hist->GetBinContent(maxBin), meanEstimate, sigmaEstimate);
838 hist->Fit("fit","RQ","Q");
840 mean = (Float_t) fit->GetParameter(1);
841 sigma = (Float_t) fit->GetParameter(2);