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) ;
97 for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
98 if (! IsValidEventSpecie(specie, list))
100 SetEventSpecie(AliRecoParam::ConvertIndex(specie)) ;
101 if ( task == AliQAv1::kRAWS ) {
102 GetRawsData(0)->SetLabelSize(0.02);
103 const Char_t *triggers[6] = {"mean", "vertex","ORA","ORC","central","semi-central"};
104 for (Int_t itr=0; itr<6; itr++) {
106 fTrEffCal[itr] = Float_t (fNumTriggersCal[itr])/Float_t (fnEventCal);
108 fTrEffPhys[itr] = Float_t (fNumTriggers[itr])/Float_t (fnEventPhys);
110 GetRawsData(197)->Fill(triggers[itr], fTrEffCal[itr]);
111 GetRawsData(197)->SetBinContent(itr+1, fTrEffCal[itr]);
112 GetRawsData(97)->Fill(triggers[itr], fTrEffPhys[itr]);
113 GetRawsData(97)->SetBinContent(itr+1, fTrEffPhys[itr]);
117 for(Int_t ik=0; ik<24; ik++)
120 effic = Float_t(feffC[ik])/Float_t(fnEventCal);
121 GetRawsData(205)->SetBinContent(ik+1,effic) ;
123 effic = Float_t(feffA[ik])/Float_t(fnEventCal);
124 GetRawsData(206)->SetBinContent(ik+1,effic );
126 effic = Float_t(feffqtc[ik])/Float_t(fnEventCal);
127 GetRawsData(207)->SetBinContent(ik+1, effic);
134 //____________________________________________________________________________
135 void AliT0QADataMakerRec::StartOfDetectorCycle()
137 //Detector specific actions at start of cycle
142 //____________________________________________________________________________
143 void AliT0QADataMakerRec::InitRaws()
145 // create Raw histograms in Raw subdir
146 const Bool_t expert = kTRUE ;
147 const Bool_t saveCorr = kTRUE ;
148 const Bool_t image = kTRUE ;
150 TString timename, ampname, qtcname, ledname;
151 TString timeCalname, ampCalname, ledCalname, qtcCalname;
152 TString qt1name, qt0name, qt1Calname, qt0Calname;
154 TH1F* fhRefPoint = new TH1F("hRefPoint","Ref Point", 10,1252170, 1252180);
155 Add2RawsList( fhRefPoint,0, expert, image, !saveCorr);
157 TH1F *fhRawCFD[24]; TH1F * fhRawLEDamp[24];
158 TH1F *fhRawQTC[24]; TH1F * fhRawLED[24];
159 TH1F *fhRawCFDcal[24]; TH1F * fhRawLEDampcal[24];
160 TH1F *fhRawCFDcalpmt[24]; TH1F * fhRawLEDcalpmt[24];
161 TH1F *fhRawCFDpmt[24]; TH1F * fhRawLEDpmt[24];
162 TH1F *fhRawQTCcal[24]; TH1F * fhRawLEDcal[24];
163 TH1F *fhRawQT0cal[24]; TH1F *fhRawQT1cal[24];
164 TH1F *fhRawQT1[24]; TH1F *fhRawQT0[24];
166 for (Int_t i=0; i<24; i++)
171 qt0name = "hRawQT0_";
172 qt1name = "hRawQT1_";
173 ampname = "hRawLEDminCFD";
180 // fhRawCFD[i] = new TH1F(timename.Data(), Form("%s;CFD [#channels];Counts", timename.Data()),10000,0,10000);
181 fhRawCFD[i] = new TH1F(timename.Data(), Form("%s;CFD [#channels];Counts", timename.Data()),10000,80000,100000);
182 Add2RawsList( fhRawCFD[i],i+1, expert, !image, !saveCorr);
183 fhRawLED[i] = new TH1F(ledname.Data(), Form("%s;LED[#channels];Counts", ledname.Data()),10000,80000,100000);
184 Add2RawsList( fhRawLED[i],i+24+1, expert, !image, !saveCorr);
185 fhRawLEDamp[i] = new TH1F(ampname.Data(), Form("%s;LED-CFD [#channels];Counts", ampname.Data()),1000,0,1000);
186 Add2RawsList( fhRawLEDamp[i],i+48+1, expert, !image, !saveCorr);
187 fhRawQTC[i] = new TH1F(qtcname.Data(), Form("%s;QTC[#channels];Counts", qtcname.Data()),10000,0,10000);
188 Add2RawsList( fhRawQTC[i],i+72+1, expert, !image, !saveCorr);
189 // fhRawQT1[i] = new TH1F(qt1name.Data(), Form("%s;QT1[#channels];Counts", qtcname.Data()),10000,80000,100000);
190 fhRawQT1[i] = new TH1F(qt1name.Data(), Form("%s;QT1[#channels];Counts", qtcname.Data()),10000,80000,100000);
191 Add2RawsList( fhRawQT1[i],270+i, expert, !image, !saveCorr);
192 fhRawQT0[i] = new TH1F(qt1name.Data(), Form("%s;QT0[#channels];Counts", qtcname.Data()),10000,80000,100000);
193 // fhRawQT0[i] = new TH1F(qt1name.Data(), Form("%s;QT0[#channels];Counts", qtcname.Data()),10000,0,10000);
194 Add2RawsList( fhRawQT0[i],270+24+i, expert, !image, !saveCorr);
196 TH1F* fhRawTrigger = new TH1F("hRawTrigger"," phys triggers;Trigger #;Counts",5,0,5);
197 Add2RawsList(fhRawTrigger ,97, !expert, !image, !saveCorr);
199 TH1F* fhRawMean = new TH1F("hRawMean","online mean signal, physics event;", 10000,80000,100000);
200 Add2RawsList( fhRawMean,98, expert, !image, !saveCorr);
201 TH1F* fhRawVertex = new TH1F("hRawVertex","online vertex signal; counts", 10000,80000,100000);
202 Add2RawsList( fhRawVertex,99, expert, !image, !saveCorr);
203 TH1F* fhRawORA = new TH1F("hRawORA","online OR A; counts", 10000,80000,100000);
204 Add2RawsList( fhRawORA,100, expert, !image, !saveCorr);
205 TH1F* fhRawORC = new TH1F("hRawORC","online OR C;counts", 10000,80000,100000);
206 Add2RawsList( fhRawORC,101, expert, !image, !saveCorr);
208 for (Int_t i=0; i<24; i++)
210 // for events with trigger CALIBRATION_EVENT
211 timeCalname ="hRawCFDcal";
212 ledCalname = "hRawLEDcal";
213 ampCalname = "hRawLEDminCFDcal";
214 qtcCalname = "hRawQTCcal";
215 qt0Calname = "hRawQT0cal";
216 qt1Calname = "hRawQT1cal";
223 // fhRawCFDcal[i] = new TH1F(timeCalname.Data(), Form("%s;Time ;Counts", timeCalname.Data()),10000,0,10000);
224 fhRawCFDcal[i] = new TH1F(timeCalname.Data(), Form("%s;Time ;Counts", timeCalname.Data()),10000,60000,80000);
225 Add2RawsList( fhRawCFDcal[i],101+i+1, expert, !image, !saveCorr);
226 // fhRawLEDcal[i] = new TH1F(ledCalname.Data(), Form("%s;Time ;Counts", ledCalname.Data()),10000,0,10000);
227 fhRawLEDcal[i] = new TH1F(ledCalname.Data(), Form("%s;Time ;Counts", ledCalname.Data()),10000,60000,80000);
228 Add2RawsList( fhRawLEDcal[i],101+i+24+1, expert, !image, !saveCorr);
229 fhRawLEDampcal[i] = new TH1F(ampCalname.Data(), Form("%s;Amplitude [ADC counts];Counts", ampCalname.Data()),1000,0,1000);
230 Add2RawsList( fhRawLEDampcal[i],101+i+48+1, expert, !image, !saveCorr);
231 fhRawQTCcal[i] = new TH1F(qtcCalname.Data(), Form("%s;Charge ;Counts",qtcCalname.Data()),10000,0,10000);
232 Add2RawsList( fhRawQTCcal[i],101+i+72+1, expert, !image, !saveCorr);
233 fhRawQT0cal[i] = new TH1F(qt0Calname.Data(), Form("%s;Charge ;Counts",qt0Calname.Data()),10000,60000,80000);
234 //fhRawQT0cal[i] = new TH1F(qt0Calname.Data(), Form("%s;Charge ;Counts",qt0Calname.Data()),10000,0,10000);
235 Add2RawsList( fhRawQT0cal[i],371+i, expert, !image, !saveCorr);
236 fhRawQT1cal[i] = new TH1F(qt1Calname.Data(), Form("%s;Charge ;Counts",qt1Calname.Data()),10000,60000,80000);
237 //fhRawQT1cal[i] = new TH1F(qt1Calname.Data(), Form("%s;Charge ;Counts",qt1Calname.Data()),10000,0,10000);
238 Add2RawsList( fhRawQT1cal[i],i+371+24, expert, !image, !saveCorr);
241 //from PMT1 (equalizing)
242 for (Int_t i=0; i<24; i++)
244 // for events with trigger CALIBRATION_EVENT
245 timeCalname ="hRawCFDcalpmt";
246 timename ="hRawCFDpmt";
251 fhRawCFDcalpmt[i] = new TH1F(timeCalname.Data(), Form("%s;Time;Counts", timeCalname.Data()),2000,-1000,1000);
252 Add2RawsList( fhRawCFDcalpmt[i],321+i , expert, !image, !saveCorr);
254 fhRawCFDpmt[i] = new TH1F(timename.Data(), Form("%s;Time;Counts", timename.Data()),2000,-1000,1000);
255 Add2RawsList( fhRawCFDcalpmt[i],220+i, expert, !image, !saveCorr);
258 TH1F* fhRawTriggerCal = new TH1F("hRawTriggerCal"," laser triggers",6,0,6);
259 Add2RawsList(fhRawTriggerCal ,197 , expert, image, !saveCorr);
261 TH1F* fhRawMeanCal = new TH1F("hRawMeanCal","online mean signal, calibration event",
263 Add2RawsList( fhRawMeanCal,198, expert, !image, !saveCorr);
264 TH1F* fhRawVertexCal = new TH1F("hRawVertexCal","online vertex signal, calibration event ", 10000,60000,80000);
265 // TH1F* fhRawVertexCal = new TH1F("hRawVertexCal","online vertex signal, calibration event ", 10000,0,10000);
266 Add2RawsList( fhRawVertexCal,199, expert, !image, !saveCorr);
267 TH1F* fhRawORAcal = new TH1F("hRawORAcal","laser OR A; counts", 10000,60000,80000);
268 // TH1F* fhRawORAcal = new TH1F("hRawORAcal","laser OR A; counts", 10000,0,10000);
269 Add2RawsList( fhRawORAcal,200, expert, !image, !saveCorr );
270 TH1F* fhRawORCcal = new TH1F("hRawORCcal","laserOR C;counts ", 10000,60000,80000);
271 // TH1F* fhRawORCcal = new TH1F("hRawORCcal","laserOR C;counts ", 10000,0,10000);
272 Add2RawsList( fhRawORCcal,201, expert, !image, !saveCorr);
274 //multiplicity trigger
275 TH1F* fhMultcal = new TH1F("hMultcal","full mulltiplicity;Multiplicity;Entries", 10000,60000,80000);
276 Add2RawsList( fhMultcal,202, expert, !image, !saveCorr );
277 TH1F* fhMultScal = new TH1F("hMultScal","full multiplicity with semi-central trigger;Multiplicity;Entries",
279 Add2RawsList( fhMultScal,203, expert, !image, !saveCorr);
280 TH1F* fhMultCcal = new TH1F("hMultCcal","full multiplicity with central trigger;Multiplicity;Entries",
282 Add2RawsList( fhMultCcal,204, expert, !image, !saveCorr);
285 TH1F* fhEffCFD = new TH1F("hEffCFD","#PMT; #CFD counts/nEvents",24, 0 ,24);
286 Add2RawsList( fhEffCFD,205, !expert, image, !saveCorr);
287 TH1F* fhEffLED = new TH1F("hEffLED","#PMT; #LED counts/nEvent",24, 0 ,24);
288 Add2RawsList( fhEffLED,206, !expert, image, !saveCorr);
289 TH1F* fhEffQTC = new TH1F("hEffQTC","#PMT; QTC efficiency%s;",24, 0 ,24);
290 Add2RawsList( fhEffQTC,207, !expert, image, !saveCorr);
292 TH2F* fhCFDcal = new TH2F("hCFDcal","CFD laser; #PMT; CFD {#channnels}",25, 0 ,25, 5000,60000,80000);
293 //TH2F* fhCFDcal = new TH2F("hCFDcal","CFD laser; #PMT; CFD {#channnels}",25, 0 ,25, 2000,0,10000);
294 fhCFDcal->SetOption("COLZ");
295 Add2RawsList( fhCFDcal,208, expert, image, !saveCorr);
296 TH2F* fhLEDcal = new TH2F("hLEDcal","LED laser; #PMT; LED [#channnels]",25, 0 ,25, 5000,60000,80000);
297 // TH2F* fhLEDcal = new TH2F("hLEDcal","LED laser; #PMT; LED [#channnels]",25, 0 ,25, 2000,0,10000);
298 fhLEDcal->SetOption("COLZ");
299 Add2RawsList( fhLEDcal,209, expert, image, !saveCorr);
301 const Char_t *triggers[6] = {"mean", "vertex","ORA","ORC","central","semi-central"};
302 for (Int_t itr=0; itr<6; itr++) {
303 GetRawsData(197)->Fill(triggers[itr], fNumTriggersCal[itr]);
304 GetRawsData(197)->SetBinContent(itr+1, fNumTriggersCal[itr]);
305 GetRawsData(97)->Fill(triggers[itr], fNumTriggers[itr]);
306 GetRawsData(97)->SetBinContent(itr+1, fNumTriggers[itr]);
310 //____________________________________________________________________________
311 void AliT0QADataMakerRec::InitDigits()
313 // create Digits histograms in Digits subdir
314 const Bool_t expert = kTRUE ;
315 const Bool_t image = kTRUE ;
317 TH2F * fhDigCFD = new TH2F("fhDigCFD", " CFD digits; #PMT; CFD digits[#channels]",25,-0.5,24.5,100,0,1000);
318 fhDigCFD->SetOption("COLZ");
319 Add2DigitsList( fhDigCFD,0, !expert, image);
320 TH2F *fhDigLEDamp = new TH2F("fhDigLEDamp", " LED-CFD digits; #PMT; LED-CFD amplitude ",25,-0.5,24.5,100,100,1000);
321 fhDigLEDamp->SetOption("COLZ");
322 Add2DigitsList( fhDigLEDamp,1, !expert, !image);
323 TH2F * fhDigQTC = new TH2F("fhDigQTC", " QTC digits; #PMT; QTC amplitude",25,-0.5,24.5,100,100,10000);
324 fhDigQTC->SetOption("COLZ");
325 Add2DigitsList( fhDigQTC,2, !expert, !image);}
327 //____________________________________________________________________________
329 void AliT0QADataMakerRec::InitRecPoints()
331 // create cluster histograms in RecPoint subdir
332 const Bool_t expert = kTRUE ;
333 const Bool_t image = kTRUE ;
335 TH2F* fhRecCFD = new TH2F("hRecCFD"," CFD time;Time [ns];Counts",24, 0 ,24,
337 fhRecCFD->SetOption("COLZ");
338 Add2RecPointsList ( fhRecCFD,0, !expert, image);
340 TH2F* fhRecAmpDiff = new TH2F("hRecAmpDiff"," LED-CFD min QTC amplitude;Amplitude [ADC counts];Counts",
341 24, 0 ,24, 200,-10,10);
342 fhRecAmpDiff->SetOption("COLZ");
343 Add2RecPointsList (fhRecAmpDiff, 1, !expert, image);
345 TH1F *fhMean = new TH1F("hMean","online - rec mean;online - rec mean[#channels];",2000, -1000, 1000);
346 Add2RecPointsList ( fhMean,2, !expert, image);
350 //____________________________________________________________________________
351 void AliT0QADataMakerRec::InitESDs()
353 //create ESDs histograms in ESDs subdir
354 const Bool_t expert = kTRUE ;
355 const Bool_t image = kTRUE ;
357 TH1F *fhESDMean = new TH1F("hESDmean"," ESD mean; mean time[%channels]",1000,0,1000);
358 Add2ESDsList(fhESDMean, 0, !expert, image) ;
359 TH1F * fhESDVertex = new TH1F("hESDvertex","ESDvertex; vertex[cm];",82,-30,30);
360 Add2ESDsList(fhESDVertex, 1, !expert, image) ;
365 //____________________________________________________________________________
366 void AliT0QADataMakerRec::MakeRaws( AliRawReader* rawReader)
370 //fills QA histos for RAW
372 Int_t refPointParam = GetRecoParam()->GetRefPoint();
375 AliT0RawReader *start = new AliT0RawReader(rawReader);
378 AliDebug(AliQAv1::GetQADebugLevel(),Form(" no raw data found!!"));
382 UInt_t type =rawReader->GetType();
383 if (type == 8){ shift=101; fnEventCal++;}
384 if (type == 7){ shift=0; fnEventPhys++;}
385 Int_t allData[110][5];
386 for (Int_t i0=0; i0<105; i0++)
388 for (Int_t j0=0; j0<5; j0++) allData[i0][j0]=0;
390 for (Int_t i=0; i<105; i++)
391 for (Int_t iHit=0; iHit<5; iHit++)
392 allData[i][iHit]= start->GetData(i,iHit);
394 if (allData[0][0] > 0 ) GetRawsData(0) -> Fill( allData[0][0]);
396 if (refPointParam < 0 ) refpoint=0;
398 refpoint = allData[0][0] - 5000;
400 for (Int_t ik = 0; ik<12; ik++){
401 for (Int_t iHt=0; iHt<5; iHt++){
403 if(allData[ik+1][iHt]>0) {
404 GetRawsData(shift+ik+1) -> Fill(allData[ik+1][iHt]-refpoint);
405 GetRawsData(shift+ik+220) -> Fill(allData[ik+1][iHt]-allData[1][0]);
408 GetRawsData(208)->Fill(ik+1, allData[ik+1][iHt]-refpoint);
412 if(allData[ik+13][iHt] > 0) {
413 GetRawsData(shift+ik+24+1)-> Fill(allData[ik+13][iHt]-refpoint);
416 GetRawsData(209)->Fill(ik+1, allData[ik+13][iHt]-refpoint);
421 if(allData[ik+13][iHt] > 0 && allData[ik+1][iHt] >0 )
422 GetRawsData(shift+ik+48+1)->
423 Fill(allData[ik+13][iHt]-allData[ik+1][iHt]);
425 if(allData[2*ik+25][iHt] > 0 || allData[2*ik+26][iHt] > 0) {
426 GetRawsData(shift+ik+72+1)->
427 Fill(allData[2*ik+26][iHt]-allData[2*ik+25][iHt]);
428 GetRawsData(shift+ik+270)->Fill(allData[2*ik+26][iHt]);
429 GetRawsData(shift+ik+24+270)->Fill(allData[2*ik+25][iHt]);
431 if(type == 8 ) feffqtc[ik]++;
435 for (Int_t ik = 12; ik<24; ik++) {
436 for (Int_t iHt=0; iHt<5; iHt++) {
437 if(allData[ik+45][iHt]>0) {
439 GetRawsData(shift+ik+1)->
440 Fill(allData[ik+45][iHt]-refpoint);
441 GetRawsData(shift+ik+220) -> Fill(allData[ik+45][iHt]-allData[57][0]);
444 GetRawsData(208)->Fill(ik+1, allData[ik+45][iHt]-refpoint);
448 if(allData[ik+57][iHt] > 0 ) {
449 GetRawsData(shift+ik+24+1)->Fill(allData[ik+57][iHt]-refpoint);
452 GetRawsData(209)->Fill(ik+1, allData[ik+57][iHt]-refpoint);
456 if(allData[2*ik+57][iHt]>0 || allData[2*ik+58][iHt]>0)
458 GetRawsData(shift+ik+72+1)-> Fill(allData[2*ik+58][iHt]-allData[2*ik+57][iHt]);
460 GetRawsData(shift+ik+270)->Fill(allData[2*ik+26][iHt]);
461 GetRawsData(shift+ik+24+270)->Fill(allData[2*ik+25][iHt]);
462 if(type == 8 ) feffqtc[ik]++;
465 if(allData[ik+57][iHt] > 0 &&allData[ik+45][iHt]>0)
466 GetRawsData(shift+ik+48+1)->
467 Fill(allData[ik+57][iHt]-allData[ik+45][iHt]);
471 Int_t trChannel[6] = {49,50,51,52,55,56};
475 for (Int_t iHt=0; iHt<6; iHt++) {
476 for (Int_t itr=0; itr<6; itr++) {
477 if(allData[trChannel[itr]][iHt]>0) {
479 GetRawsData(98+itr)->Fill(allData[trChannel[itr]][iHt]);
487 for (Int_t iHt=0; iHt<5; iHt++) {
488 for (Int_t itr=0; itr<6; itr++) {
489 if(allData[trChannel[itr]][iHt]>0)
492 GetRawsData(198+itr)->Fill(allData[trChannel[itr]][iHt]);
493 fNumTriggersCal[itr]++;
496 if(allData[53][iHt]>0 && allData[54][iHt]>0) {
497 GetRawsData(204)->Fill(allData[53][iHt]-allData[54][iHt]);
498 if(allData[55][iHt]>0) GetRawsData(202)->Fill(allData[53][iHt]-allData[54][iHt]);
499 if(allData[56][iHt]>0) GetRawsData(203)->Fill(allData[53][iHt]-allData[54][iHt]);
510 //____________________________________________________________________________
511 void AliT0QADataMakerRec::MakeDigits( TTree *digitsTree)
513 //fills QA histos for Digits
515 TArrayI *digCFD = new TArrayI(24);
516 TArrayI *digLED = new TArrayI(24);
517 TArrayI *digQT0 = new TArrayI(24);
518 TArrayI *digQT1 = new TArrayI(24);
521 TBranch *brDigits=digitsTree->GetBranch("T0");
522 AliT0digit *fDigits = new AliT0digit() ;
524 brDigits->SetAddress(&fDigits);
526 AliError(Form("EXEC Branch T0 digits not found"));
529 digitsTree->GetEvent(0);
530 digitsTree->GetEntry(0);
531 brDigits->GetEntry(0);
532 fDigits->GetTimeCFD(*digCFD);
533 fDigits->GetTimeLED(*digLED);
534 fDigits->GetQT0(*digQT0);
535 fDigits->GetQT1(*digQT1);
536 refpoint = fDigits->RefPoint();
537 for (Int_t i=0; i<24; i++)
539 if (digCFD->At(i)>0) {
540 Int_t cfd=digCFD->At(i)- refpoint;
541 GetDigitsData(0) ->Fill(i,cfd);
542 GetDigitsData(1) -> Fill(i, (digLED->At(i) - digCFD->At(i)));
543 GetDigitsData(2) -> Fill(i, (digQT1->At(i) - digQT0->At(i)));
554 //____________________________________________________________________________
555 void AliT0QADataMakerRec::MakeRecPoints(TTree * clustersTree)
557 //fills QA histos for clusters
559 AliT0RecPoint* frecpoints= new AliT0RecPoint ();
561 AliError(":MakeRecPoints >> no recpoints found");
564 TBranch *brRec =clustersTree ->GetBranch("T0");
566 brRec->SetAddress(&frecpoints);
568 AliError(Form("EXEC Branch T0 rec not found "));
574 for ( Int_t i=0; i<24; i++) {
576 GetRecPointsData(0) -> Fill(i, frecpoints -> GetTime(i) - frecpoints -> GetTime(0));
578 GetRecPointsData(0) -> Fill(i, frecpoints -> GetTime(i) - frecpoints -> GetTime(12));
579 GetRecPointsData(1) -> Fill( i, frecpoints -> GetAmp(i) - frecpoints->AmpLED(i));
581 Double_t mmm=frecpoints->GetOnlineMean()- frecpoints->GetMeanTime();
582 GetRecPointsData(2) ->Fill(mmm);
586 //____________________________________________________________________________
587 void AliT0QADataMakerRec::MakeESDs(AliESDEvent * esd)
589 //fills QA histos for ESD
591 GetESDsData(0) -> Fill(esd->GetT0());
592 GetESDsData(1)-> Fill(esd->GetT0zVertex());