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>
30 // --- Standard library ---
32 // --- AliRoot header files ---
33 #include "AliESDEvent.h"
35 #include "AliT0digit.h"
37 #include "AliT0RecPoint.h"
38 #include "AliT0QADataMakerRec.h"
39 #include "AliQAChecker.h"
40 #include "AliT0RawReader.h"
41 #include "AliT0RecoParam.h"
43 #include "Riostream.h"
44 ClassImp(AliT0QADataMakerRec)
46 //____________________________________________________________________________
47 AliT0QADataMakerRec::AliT0QADataMakerRec() :
48 AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kT0),
49 "T0 Quality Assurance Data Maker"),
55 for (Int_t i=0; i<6; i++) {
62 for (Int_t i=0; i<24; i++)
70 //____________________________________________________________________________
71 AliT0QADataMakerRec::AliT0QADataMakerRec(const AliT0QADataMakerRec& qadm) :
79 SetName((const char*)qadm.GetName()) ;
80 SetTitle((const char*)qadm.GetTitle());
83 //__________________________________________________________________
84 AliT0QADataMakerRec& AliT0QADataMakerRec::operator = (const AliT0QADataMakerRec& qadm )
87 this->~AliT0QADataMakerRec();
88 new(this) AliT0QADataMakerRec(qadm);
91 //____________________________________________________________________________
92 void AliT0QADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
94 //Detector specific actions at end of cycle
96 AliQAChecker::Instance()->Run(AliQAv1::kT0, task, list) ;
98 for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
99 if (! IsValidEventSpecie(specie, list))
101 SetEventSpecie(AliRecoParam::ConvertIndex(specie)) ;
102 if ( task == AliQAv1::kRAWS ) {
103 GetRawsData(0)->SetLabelSize(0.02);
104 const Char_t *triggers[6] = {"mean", "vertex","ORA","ORC","central","semi-central"};
105 for (Int_t itr=0; itr<6; itr++) {
107 fTrEffCal[itr] = Float_t (fNumTriggersCal[itr])/Float_t (fnEventCal);
109 fTrEffPhys[itr] = Float_t (fNumTriggers[itr])/Float_t (fnEventPhys);
111 GetRawsData(420)->Fill(triggers[itr], fTrEffCal[itr]);
112 GetRawsData(420)->SetBinContent(itr+1, fTrEffCal[itr]);
113 GetRawsData(97)->Fill(triggers[itr], fTrEffPhys[itr]);
114 GetRawsData(97)->SetBinContent(itr+1, fTrEffPhys[itr]);
117 for(Int_t ik=0; ik<24; ik++)
120 if ( fnEventCal>0) effic = Float_t(feffC[ik])/Float_t(fnEventCal);
121 GetRawsData(205)->SetBinContent(ik+1,effic) ;
123 if ( fnEventCal>0) effic = Float_t(feffA[ik])/Float_t(fnEventCal);
124 GetRawsData(206)->SetBinContent(ik+1,effic );
126 if ( fnEventCal>0) effic = Float_t(feffqtc[ik])/Float_t(fnEventCal);
127 GetRawsData(207)->SetBinContent(ik+1, effic);
136 //____________________________________________________________________________
137 void AliT0QADataMakerRec::StartOfDetectorCycle()
139 //Detector specific actions at start of cycle
144 //____________________________________________________________________________
145 void AliT0QADataMakerRec::InitRaws()
147 // create Raw histograms in Raw subdir
148 const Bool_t expert = kTRUE ;
149 const Bool_t saveCorr = kTRUE ;
150 const Bool_t image = kTRUE ;
155 for (Int_t i=0; i<500; i++){
156 // low[i] = GetRecoParam()->GetLow(i);
157 /// high[i] = GetRecoParam()->GetHigh(i);
163 TString timename, ampname, qtcname, ledname;
164 TString timeCalname, ampCalname, ledCalname, qtcCalname;
165 TString qt1name, qt0name, qt1Calname, qt0Calname;
168 TH1F* fhRefPoint = new TH1F("hRefPoint","Ref Point",high[0]-low[0],low[0],high[0]);
169 Add2RawsList( fhRefPoint,0, expert, !image, !saveCorr);
171 TH1F *fhRawCFD[24]; TH1F * fhRawLEDamp[24];
172 TH1F *fhRawQTC[24]; TH1F * fhRawLED[24];
173 TH1F *fhRawCFDcal[24]; TH1F * fhRawLEDampcal[24];
174 TH1F *fhRawCFDcalpmt[24];
175 TH1F *fhRawCFDpmt[24];
176 TH1F *fhRawQTCcal[24]; TH1F * fhRawLEDcal[24];
177 TH1F *fhRawQT0cal[24]; TH1F *fhRawQT1cal[24];
178 TH1F *fhRawQT1[24]; TH1F *fhRawQT0[24];
179 TH1F* fhRawNhits[24];
181 for (Int_t i=0; i<24; i++)
186 qt0name = "hRawQT0_";
187 qt1name = "hRawQT1_";
188 ampname = "hRawLEDminCFD";
197 fhRawCFD[i] = new TH1F(timename.Data(), Form("%s;CFD [#channels];Counts", timename.Data()),high[i+1]-low[i+1],low[i+1],high[i+1]);
198 Add2RawsList( fhRawCFD[i],i+1, expert, !image, !saveCorr);
199 fhRawLED[i] = new TH1F(ledname.Data(), Form("%s;LED[#channels];Counts", ledname.Data()),high[i+24+1]-low[i+24+1],low[i+24+1],high[i]);
200 Add2RawsList( fhRawLED[i],i+24+1, expert, !image, !saveCorr);
201 fhRawLEDamp[i] = new TH1F(ampname.Data(), Form("%s;LED-CFD [#channels];Counts", ampname.Data()),1000,0,1000);
202 Add2RawsList( fhRawLEDamp[i],i+48+1, expert, !image, !saveCorr);
203 fhRawQTC[i] = new TH1F(qtcname.Data(), Form("%s;QTC[#channels];Counts", qtcname.Data()),10000,0,10000);
204 Add2RawsList( fhRawQTC[i],i+72+1, expert, !image, !saveCorr);
205 fhRawQT1[i] = new TH1F(qt1name.Data(), Form("%s;QT1[#channels];Counts", qtcname.Data()),high[270+i]-low[270+i],low[270+i],high[270+i]);
206 Add2RawsList( fhRawQT1[i],270+i, expert, !image, !saveCorr);
207 fhRawQT0[i] = new TH1F(qt0name.Data(), Form("%s;QT0[#channels];Counts", qtcname.Data()),high[270+24+i]-low[270+24+i],low[270+24+i],high[270+24+i]);
208 Add2RawsList( fhRawQT0[i],270+24+i, expert, !image, !saveCorr);
210 fhRawNhits[i] = new TH1F(nhits.Data(), Form("%s;#Hits;Events", nhits.Data()),10, 0, 10);
211 Add2RawsList( fhRawNhits[i],244+i, expert, !image, !saveCorr);
214 TH1F* fhRawTrigger = new TH1F("hRawTrigger"," phys triggers;Trigger #;Counts",5,0,5);
215 Add2RawsList(fhRawTrigger ,97, expert, !image, !saveCorr);
217 TH1F* fhRawMean = new TH1F("hRawMean","online mean signal, physics event;",high[98]-low[98],low[98],high[98]);
218 Add2RawsList( fhRawMean,98, expert, !image, !saveCorr);
220 TH1F* fhRawVertex = new TH1F("hRawVertex","online vertex signal; counts",high[100]-low[99],low[99],high[99]);
221 Add2RawsList( fhRawVertex,99, expert, !image, !saveCorr);
223 TH1F* fhRawORA = new TH1F("hRawORA","online OR A; counts", high[100]-low[100],low[100],high[100]);
224 Add2RawsList( fhRawORA,100, expert, !image, !saveCorr);
225 TH1F* fhRawORC = new TH1F("hRawORC","online OR C;counts", high[101]-low[101],low[101],high[101]);
226 Add2RawsList( fhRawORC,101, expert, !image, !saveCorr);
229 for (Int_t i=0; i<24; i++)
231 // for events with trigger CALIBRATION_EVENT
232 timeCalname ="hRawCFDcal";
233 ledCalname = "hRawLEDcal";
234 ampCalname = "hRawLEDminCFDcal";
235 qtcCalname = "hRawQTCcal";
236 qt0Calname = "hRawQT0cal";
237 qt1Calname = "hRawQT1cal";
245 fhRawCFDcal[i] = new TH1F(timeCalname.Data(), Form("%s;Time ;Counts", timeCalname.Data()),high[101+i+1]-low[101+i+1],low[101+i+1],high[101+i+1]);
246 Add2RawsList( fhRawCFDcal[i],101+i+1, expert, !image, !saveCorr);
248 fhRawLEDcal[i] = new TH1F(ledCalname.Data(), Form("%s;Time ;Counts", ledCalname.Data()),high[101+i+24+1]-low[101+i+24+1],low[101+i+24+1],high[101+i+24+1]);
249 Add2RawsList( fhRawLEDcal[i],101+i+24+1, expert, !image, !saveCorr);
251 fhRawLEDampcal[i] = new TH1F(ampCalname.Data(), Form("%s;Amplitude [ADC counts];Counts", ampCalname.Data()),1000,0,1000);
252 Add2RawsList( fhRawLEDampcal[i],101+i+48+1, expert, !image, !saveCorr);
254 fhRawQTCcal[i] = new TH1F(qtcCalname.Data(), Form("%s;Charge ;Counts",qtcCalname.Data()),10000,0,10000);
255 Add2RawsList( fhRawQTCcal[i],101+i+72+1, expert, !image, !saveCorr);
257 fhRawQT0cal[i] = new TH1F(qt0Calname.Data(), Form("%s;Charge ;Counts",qt0Calname.Data()),high[371+i]-low[371+i],low[371+i],high[371+i]);
258 Add2RawsList( fhRawQT0cal[i],371+i, expert, !image, !saveCorr);
260 fhRawQT1cal[i] = new TH1F(qt1Calname.Data(), Form("%s;Charge ;Counts",qt1Calname.Data()),high[i+371+24]-low[i+371+24],low[i+371+24],high[i+371+24]);
261 Add2RawsList( fhRawQT1cal[i],i+371+24, expert, !image, !saveCorr);
265 //from PMT1 (equalizing)
266 for (Int_t i=0; i<24; i++)
268 // for events with trigger CALIBRATION_EVENT
269 timeCalname ="hRawCFDcalpmt";
270 timename ="hRawCFDpmt";
275 fhRawCFDcalpmt[i] = new TH1F(timeCalname.Data(), Form("%s;Time;Counts", timeCalname.Data()),2000,-1000,1000);
276 Add2RawsList( fhRawCFDcalpmt[i],321+i , expert, !image, !saveCorr);
278 fhRawCFDpmt[i] = new TH1F(timename.Data(), Form("%s;Time;Counts", timename.Data()),2000,-1000,1000);
279 Add2RawsList( fhRawCFDpmt[i],220+i, expert, !image, !saveCorr);
282 TH1F* fhRawTriggerCal = new TH1F("hRawTriggerCal"," laser triggers",6,0,6);
283 Add2RawsList(fhRawTriggerCal ,420 , !expert, image, !saveCorr);
285 TH1F* fhRawMeanCal = new TH1F("hRawMeanCal","online mean signal, calibration event",high[198]-low[198],low[198],high[198]);
286 Add2RawsList( fhRawMeanCal,198, expert, !image, !saveCorr);
288 TH1F* fhRawVertexCal = new TH1F("hRawVertexCal","online vertex signal, calibration event ", high[199]-low[199],low[199],high[199] );
289 Add2RawsList( fhRawVertexCal,199, expert, !image, !saveCorr);
292 TH1F* fhRawORAcal = new TH1F("hRawORAcal","laser OR A; counts", high[200]-low[200],low[200],high[200]);
293 Add2RawsList( fhRawORAcal,200, expert, !image, !saveCorr );
296 TH1F* fhRawORCcal = new TH1F("hRawORCcal","laserOR C;counts ", high[201]-low[201],low[201],high[201]);
297 Add2RawsList( fhRawORCcal,201, expert, !image, !saveCorr);
300 //multiplicity trigger
301 TH1F* fhMultcal = new TH1F("hMultcal","full mulltiplicity;Multiplicity;Entries", high[202]-low[202],low[202],high[202]);
302 Add2RawsList( fhMultcal,202, expert, !image, !saveCorr );
303 TH1F* fhMultScal = new TH1F("hMultSemical","full multiplicity with semi-central trigger;Multiplicity;Entries",
304 high[203]-low[203],low[203],high[203] );
305 Add2RawsList( fhMultScal,203, expert, !image, !saveCorr);
306 TH1F* fhMultCcal = new TH1F("hMultCentrcal","full multiplicity with central trigger;Multiplicity;Entries",
307 high[204]-low[204],low[204],high[204]);
308 Add2RawsList( fhMultCcal,204, expert, !image, !saveCorr);
310 TH1F* fhMult = new TH1F("hMult","full mulltiplicity;Multiplicity;Entries", high[216]-low[216],low[216],high[216]);
311 Add2RawsList( fhMult,216, expert, !image, !saveCorr );
312 TH1F* fhMultS = new TH1F("hMultSemi","full multiplicity with semi-central trigger;Multiplicity;Entries",
313 high[217]-low[217],low[217],high[217] );
314 Add2RawsList( fhMultS,217, expert, !image, !saveCorr);
315 TH1F* fhMultC = new TH1F("hMultCentr","full multiplicity with central trigger;Multiplicity;Entries",
316 high[218]-low[218],low[218],high[218]);
317 Add2RawsList( fhMultC,218, expert, !image, !saveCorr);
319 TH1F* fhEffCFD = new TH1F("hEffCFD","#PMT; #CFD counts/nEvents",24, 0 ,24);
320 Add2RawsList( fhEffCFD,205, !expert, image, !saveCorr);
322 TH1F* fhEffLED = new TH1F("hEffLED","#PMT; #LED counts/nEvent",24, 0 ,24);
323 Add2RawsList( fhEffLED,206, !expert, image, !saveCorr);
325 TH1F* fhEffQTC = new TH1F("hEffQTC","#PMT; QTC efficiency%s;",24, 0 ,24);
326 Add2RawsList( fhEffQTC,207, !expert, image, !saveCorr);
329 TH2F* fhCFDcal = new TH2F("hCFDcal","CFD laser; #PMT; CFD {#channnels}",25, 0 ,25,high[208]-low[208],low[208],high[208]);
330 fhCFDcal->SetOption("COLZ");
331 Add2RawsList( fhCFDcal,208, !expert, image, !saveCorr);
334 TH2F* fhLEDcal = new TH2F("hLEDcal","LED laser; #PMT; LED [#channnels]",25, 0 ,25, high[209]-low[209],low[209],high[209]);
335 fhLEDcal->SetOption("COLZ");
336 Add2RawsList( fhLEDcal,209, !expert, image, !saveCorr);
339 TH1F* fhNumPMTA= new TH1F("hNumPMTA","number of PMT hitted per event",12, 0 ,12);
340 Add2RawsList(fhNumPMTA ,211, expert, !image, !saveCorr);
342 TH1F* fhNumPMTC= new TH1F("hNumPMTC","number of PMT hitted per event",12, 0 ,12);
343 Add2RawsList(fhNumPMTC ,212, expert, !image, !saveCorr);
345 TH1F* fhHitsOrA= new TH1F("fhHitsOrA","T0_OR A hit multiplicitie",10, 0 ,10);
346 Add2RawsList( fhHitsOrA,213, expert, !image, !saveCorr);
348 TH1F* fhHitsOrC= new TH1F("fhHitsOrC","T0_OR C hit multiplicitie",10, 0 ,10);
349 Add2RawsList(fhHitsOrC ,214, expert, !image, !saveCorr);
351 TH1F* fhOrCminOrA= new TH1F("fhOrCminOrA","T0_OR C - T0_OR A", high[215]-low[215],low[215],high[215]);
352 Add2RawsList( fhOrCminOrA,215, expert, !image, !saveCorr);
354 const Char_t *triggers[6] = {"mean", "vertex","ORA","ORC","central","semi-central"};
355 for (Int_t itr=0; itr<6; itr++) {
356 GetRawsData(420)->Fill(triggers[itr], fNumTriggersCal[itr]);
357 GetRawsData(420)->SetBinContent(itr+1, fNumTriggersCal[itr]);
358 GetRawsData(97)->Fill(triggers[itr], fNumTriggers[itr]);
359 GetRawsData(97)->SetBinContent(itr+1, fNumTriggers[itr]);
363 //____________________________________________________________________________
364 void AliT0QADataMakerRec::InitDigits()
366 // create Digits histograms in Digits subdir
367 const Bool_t expert = kTRUE ;
368 const Bool_t image = kTRUE ;
370 TH2F * fhDigCFD = new TH2F("fhDigCFD", " CFD digits; #PMT; CFD digits[#channels]",25,-0.5,24.5,100,0,1000);
371 fhDigCFD->SetOption("COLZ");
372 Add2DigitsList( fhDigCFD,0, !expert, image);
373 TH2F *fhDigLEDamp = new TH2F("fhDigLEDamp", " LED-CFD digits; #PMT; LED-CFD amplitude ",25,-0.5,24.5,100,100,1000);
374 fhDigLEDamp->SetOption("COLZ");
375 Add2DigitsList( fhDigLEDamp,1, !expert, !image);
376 TH2F * fhDigQTC = new TH2F("fhDigQTC", " QTC digits; #PMT; QTC amplitude",25,-0.5,24.5,100,100,10000);
377 fhDigQTC->SetOption("COLZ");
378 Add2DigitsList( fhDigQTC,2, !expert, !image);}
380 //____________________________________________________________________________
382 void AliT0QADataMakerRec::InitRecPoints()
384 // create cluster histograms in RecPoint subdir
385 const Bool_t expert = kTRUE ;
386 const Bool_t image = kTRUE ;
388 TH2F* fhRecCFD = new TH2F("hRecCFD"," CFD time;Time [ns];Counts",24, 0 ,24,
390 fhRecCFD->SetOption("COLZ");
391 Add2RecPointsList ( fhRecCFD,0, !expert, image);
393 TH2F* fhRecAmpDiff = new TH2F("hRecAmpDiff"," LED-CFD min QTC amplitude;Amplitude [ADC counts];Counts",
394 24, 0 ,24, 200,-10,10);
395 fhRecAmpDiff->SetOption("COLZ");
396 Add2RecPointsList (fhRecAmpDiff, 1, !expert, image);
398 TH1F *fhMean = new TH1F("hMean","online - rec mean;online - rec mean[#channels];",2000, -1000, 1000);
399 Add2RecPointsList ( fhMean,2, !expert, image);
403 //____________________________________________________________________________
404 void AliT0QADataMakerRec::InitESDs()
406 //create ESDs histograms in ESDs subdir
407 const Bool_t expert = kTRUE ;
408 const Bool_t image = kTRUE ;
410 TH1F *fhESDMean = new TH1F("hESDmean"," ESD mean; mean time[%channels]",1000,0,1000);
411 Add2ESDsList(fhESDMean, 0, !expert, image) ;
412 TH1F * fhESDVertex = new TH1F("hESDvertex","ESDvertex; vertex[cm];",82,-30,30);
413 Add2ESDsList(fhESDVertex, 1, !expert, image) ;
418 //____________________________________________________________________________
419 void AliT0QADataMakerRec::MakeRaws( AliRawReader* rawReader)
423 //fills QA histos for RAW
425 // Int_t refPointParam = GetRecoParam()->GetRefPoint();
427 Int_t refPointParam = 0;
429 AliT0RawReader *start = new AliT0RawReader(rawReader);
432 AliDebug(AliQAv1::GetQADebugLevel(),Form(" no raw data found!!"));
436 UInt_t type =rawReader->GetType();
437 if (type == 8){ shift=101; fnEventCal++;}
438 if (type == 7){ shift=0; fnEventPhys++;}
439 Int_t allData[110][5];
440 for (Int_t i0=0; i0<105; i0++)
442 for (Int_t j0=0; j0<5; j0++) allData[i0][j0]=0;
444 for (Int_t i=0; i<105; i++)
445 for (Int_t iHit=0; iHit<5; iHit++)
446 allData[i][iHit]= start->GetData(i,iHit);
448 if (allData[0][0] > 0 ) GetRawsData(0) -> Fill( allData[0][0]);
450 refpoint = allData[refPointParam][0];
451 if (refPointParam < 0 ) refpoint=0;
452 if (refPointParam == 0 ) refpoint = allData[0][0] - 5000;
455 for (Int_t ik = 0; ik<12; ik++)
457 // for (Int_t iHt=0; iHt<1; iHt++){
460 if(allData[ik+1][0]>0 && type == 7 ) numPmtC++;
461 for (Int_t iHt=0; iHt<5; iHt++){
463 if(allData[ik+1][iHt]>0) {
464 GetRawsData(shift+ik+1) -> Fill(allData[ik+1][iHt]-refpoint);
466 GetRawsData(shift+ik+220) -> Fill(allData[ik+1][iHt]-allData[1][0]);
467 // cout<<"C cfd "<<ik<<" "<<iHt<<" "<<allData[ik+1][iHt]<<" -RF "<<refpoint<<" "<<allData[ik+1][iHt]-refpoint<<endl;
470 GetRawsData(208)->Fill(ik+1, allData[ik+1][iHt]-refpoint);
472 if(type == 7 ) nhitsPMT++;
476 if(allData[ik+13][iHt] > 0) {
477 GetRawsData(shift+ik+24+1)-> Fill(allData[ik+13][iHt]-refpoint);
478 // cout<<"C led "<<ik<<" "<<iHt<<" "<<allData[ik+1][iHt]<<" -RF "<<refpoint<<" "<<allData[ik+13][iHt]-refpoint<<endl;
481 GetRawsData(209)->Fill(ik+1, allData[ik+13][iHt]-refpoint);
486 if(allData[ik+13][iHt] > 0 && allData[ik+1][iHt] >0 )
487 GetRawsData(shift+ik+48+1)->
488 Fill(allData[ik+13][iHt]-allData[ik+1][iHt]);
491 if(allData[2*ik+25][iHt] > 0 || allData[2*ik+26][iHt] > 0) {
492 GetRawsData(shift+ik+72+1)->
493 Fill(allData[2*ik+25][iHt]-allData[2*ik+26][iHt]);
494 GetRawsData(shift+ik+270)->Fill(allData[2*ik+26][iHt] - refpoint);
495 GetRawsData(shift+ik+24+270)->Fill(allData[2*ik+25][iHt] - refpoint);
497 if(type == 8 ) feffqtc[ik]++;
500 if(type == 7 ) GetRawsData(ik+244)->Fill(nhitsPMT);
502 if(type == 7 ) GetRawsData(212)->Fill(numPmtC);
504 for (Int_t ik = 12; ik<24; ik++) {
505 // for (Int_t iHt=0; iHt<1; iHt++) {
506 if(allData[ik+45][0]>0 && type == 7 ) numPmtA++;
508 for (Int_t iHt=0; iHt<5; iHt++) {
509 if(allData[ik+45][iHt]>0) {
511 GetRawsData(shift+ik+1)->
512 Fill(allData[ik+45][iHt]-refpoint);
514 GetRawsData(shift+ik+220) -> Fill(allData[ik+45][iHt]-allData[57][0]);
516 // cout<<"A cfd "<<ik<<" "<<iHt<<" "<<allData[ik+1][iHt]<<" -RF "<<refpoint<<" "<<allData[ik+1][iHt]-refpoint<<endl;
519 GetRawsData(208)->Fill(ik+1, allData[ik+45][iHt]-refpoint);
521 if(type == 7 ) nhitsPMT++;
526 if(allData[ik+57][iHt] > 0 ) {
527 GetRawsData(shift+ik+24+1)->Fill(allData[ik+57][iHt]-refpoint);
528 // cout<<"C led "<<ik<<" "<<iHt<<" "<<allData[ik+1][iHt]<<" -RF "<<refpoint<<" "<<allData[ik+13][iHt]-refpoint<<endl;
531 GetRawsData(209)->Fill(ik+1, allData[ik+57][iHt]-refpoint);
535 if(allData[2*ik+57][iHt]>0 || allData[2*ik+58][iHt]>0)
537 GetRawsData(shift+ik+72+1)-> Fill(allData[2*ik+57][iHt]-allData[2*ik+58][iHt]);
538 GetRawsData(shift+ik+270)->Fill(allData[2*ik+58][iHt]-refpoint);
539 GetRawsData(shift+ik+24+270)->Fill(allData[2*ik+57][iHt]-refpoint);
540 if(type == 8 ) feffqtc[ik]++;
543 if(allData[ik+57][iHt] > 0 &&allData[ik+45][iHt]>0)
544 GetRawsData(shift+ik+48+1)->
545 Fill(allData[ik+57][iHt]-allData[ik+45][iHt]);
547 if(type == 7 ) GetRawsData(ik+244)->Fill(nhitsPMT);
550 if(type == 7 ) GetRawsData(211)->Fill(numPmtA);
552 Int_t trChannel[6] = {49,50,51,52,55,56};
558 for (Int_t iHt=0; iHt<5; iHt++) {
559 for (Int_t itr=0; itr<6; itr++) {
560 if (allData[trChannel[itr]][iHt] >0) {
562 GetRawsData(98+itr)->Fill(allData[trChannel[itr]][iHt]-refpoint);
566 if(allData[51][iHt] >0) nhitsOrA++;
567 if(allData[52][iHt] >0) nhitsOrC++;
569 if(allData[53][iHt]>0 && allData[54][iHt]>0) {
570 GetRawsData(216)->Fill(allData[53][iHt]-allData[54][iHt]);
571 if(allData[55][iHt]>0) GetRawsData(216)->Fill(allData[53][iHt]-allData[54][iHt]);
572 if(allData[56][iHt]>0) GetRawsData(218)->Fill(allData[53][iHt]-allData[54][iHt]);
574 if(allData[51][iHt]>0 && allData[52][iHt]>0)
575 GetRawsData(215)->Fill(allData[52][iHt]-allData[51][iHt]);
577 GetRawsData(213)->Fill(nhitsOrA);
578 GetRawsData(214)->Fill(nhitsOrC);
584 for (Int_t iHt=0; iHt<5; iHt++) {
585 for (Int_t itr=0; itr<6; itr++) {
586 if(allData[trChannel[itr]][iHt]>0)
589 GetRawsData(198+itr)->
590 Fill(allData[trChannel[itr]][iHt]-refpoint);
591 fNumTriggersCal[itr]++;
594 if(allData[53][iHt]>0 && allData[54][iHt]>0) {
595 GetRawsData(202)->Fill(allData[53][iHt]-allData[54][iHt]);
596 if(allData[55][iHt]>0) GetRawsData(202)->Fill(allData[53][iHt]-allData[54][iHt]);
597 if(allData[56][iHt]>0) GetRawsData(204)->Fill(allData[53][iHt]-allData[54][iHt]);
608 //____________________________________________________________________________
609 void AliT0QADataMakerRec::MakeDigits( TTree *digitsTree)
611 //fills QA histos for Digits
613 TArrayI *digCFD = new TArrayI(24);
614 TArrayI *digLED = new TArrayI(24);
615 TArrayI *digQT0 = new TArrayI(24);
616 TArrayI *digQT1 = new TArrayI(24);
619 TBranch *brDigits=digitsTree->GetBranch("T0");
620 AliT0digit *fDigits = new AliT0digit() ;
622 brDigits->SetAddress(&fDigits);
624 AliError(Form("EXEC Branch T0 digits not found"));
627 digitsTree->GetEvent(0);
628 digitsTree->GetEntry(0);
629 brDigits->GetEntry(0);
630 fDigits->GetTimeCFD(*digCFD);
631 fDigits->GetTimeLED(*digLED);
632 fDigits->GetQT0(*digQT0);
633 fDigits->GetQT1(*digQT1);
634 refpoint = fDigits->RefPoint();
635 for (Int_t i=0; i<24; i++)
637 if (digCFD->At(i)>0) {
638 Int_t cfd=digCFD->At(i)- refpoint;
639 GetDigitsData(0) ->Fill(i,cfd);
640 GetDigitsData(1) -> Fill(i, (digLED->At(i) - digCFD->At(i)));
641 GetDigitsData(2) -> Fill(i, (digQT1->At(i) - digQT0->At(i)));
652 //____________________________________________________________________________
653 void AliT0QADataMakerRec::MakeRecPoints(TTree * clustersTree)
655 //fills QA histos for clusters
657 AliT0RecPoint* frecpoints= new AliT0RecPoint ();
659 AliError(":MakeRecPoints >> no recpoints found");
662 TBranch *brRec =clustersTree ->GetBranch("T0");
664 brRec->SetAddress(&frecpoints);
666 AliError(Form("EXEC Branch T0 rec not found "));
672 for ( Int_t i=0; i<24; i++) {
674 GetRecPointsData(0) -> Fill(i, frecpoints -> GetTime(i) - frecpoints -> GetTime(0));
676 GetRecPointsData(0) -> Fill(i, frecpoints -> GetTime(i) - frecpoints -> GetTime(12));
677 GetRecPointsData(1) -> Fill( i, frecpoints -> GetAmp(i) - frecpoints->AmpLED(i));
679 Double_t mmm=frecpoints->GetOnlineMean()- frecpoints->GetMeanTime();
680 GetRecPointsData(2) ->Fill(mmm);
684 //____________________________________________________________________________
685 void AliT0QADataMakerRec::MakeESDs(AliESDEvent * esd)
687 //fills QA histos for ESD
689 GetESDsData(0) -> Fill(esd->GetT0());
690 GetESDsData(1)-> Fill(esd->GetT0zVertex());