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"
42 #include "Riostream.h"
43 ClassImp(AliT0QADataMakerRec)
45 //____________________________________________________________________________
46 AliT0QADataMakerRec::AliT0QADataMakerRec() :
47 AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kT0),
48 "T0 Quality Assurance Data Maker"),
53 for (Int_t i=0; i<6; i++) {
57 for (Int_t i=0; i<24; i++)
65 //____________________________________________________________________________
66 AliT0QADataMakerRec::AliT0QADataMakerRec(const AliT0QADataMakerRec& qadm) :
72 SetName((const char*)qadm.GetName()) ;
73 SetTitle((const char*)qadm.GetTitle());
76 //__________________________________________________________________
77 AliT0QADataMakerRec& AliT0QADataMakerRec::operator = (const AliT0QADataMakerRec& qadm )
80 this->~AliT0QADataMakerRec();
81 new(this) AliT0QADataMakerRec(qadm);
84 //____________________________________________________________________________
85 void AliT0QADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
87 //Detector specific actions at end of cycle
89 AliQAChecker::Instance()->Run(AliQAv1::kT0, task, list) ;
90 for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
91 if (! IsValidEventSpecie(specie, list))
93 SetEventSpecie(AliRecoParam::ConvertIndex(specie)) ;
94 if ( task == AliQAv1::kRAWS ) {
95 const Char_t *triggers[6] = {"mean", "vertex","ORA","ORC","central","semi-central"};
96 for (Int_t itr=0; itr<6; itr++) {
97 GetRawsData(197)->Fill(triggers[itr], fNumTriggersCal[itr]);
98 GetRawsData(197)->SetBinContent(itr+1, fNumTriggersCal[itr]);
104 //____________________________________________________________________________
105 void AliT0QADataMakerRec::StartOfDetectorCycle()
107 //Detector specific actions at start of cycle
112 //____________________________________________________________________________
113 void AliT0QADataMakerRec::InitRaws()
115 // create Raw histograms in Raw subdir
116 const Bool_t expert = kTRUE ;
117 const Bool_t saveCorr = kTRUE ;
118 const Bool_t image = kTRUE ;
120 TString timename, ampname, qtcname, ledname;
121 TString timeCalname, ampCalname, ledCalname, qtcCalname;
123 TH1F* fhRefPoint = new TH1F("hRefPoint","Ref Point", 10,1252170, 1252180);
124 Add2RawsList( fhRefPoint,0, !expert, image, !saveCorr);
126 TH1F *fhRawCFD[24]; TH1F * fhRawLEDamp[24];
127 TH1F *fhRawQTC[24]; TH1F * fhRawLED[24];
128 TH1F *fhRawCFDcal[24]; TH1F * fhRawLEDampcal[24];
129 TH1F *fhRawQTCcal[24]; TH1F * fhRawLEDcal[24];
131 for (Int_t i=0; i<24; i++)
136 ampname = "hRawLEDminCFD";
141 fhRawCFD[i] = new TH1F(timename.Data(), Form("%s;CFD [#channels];Counts", timename.Data()),10000,0,10000);
142 Add2RawsList( fhRawCFD[i],i+1, expert, !image, !saveCorr);
143 fhRawLED[i] = new TH1F(ledname.Data(), Form("%s;LED[#channels];Counts", ledname.Data()),10000,0,10000);
144 Add2RawsList( fhRawLED[i],i+24+1, expert, !image, !saveCorr);
145 fhRawLEDamp[i] = new TH1F(ampname.Data(), Form("%s;LED-CFD [#channels];Counts", ampname.Data()),10000,0,10000);
146 Add2RawsList( fhRawLEDamp[i],i+48+1, expert, !image, !saveCorr);
147 fhRawQTC[i] = new TH1F(qtcname.Data(), Form("%s;QTC[#channels];Counts", qtcname.Data()),700,0,7000);
148 Add2RawsList( fhRawQTC[i],i+72+1, expert, image, !saveCorr);
150 TH1F* fhRawTrigger = new TH1F("hRawTrigger"," phys triggers;Trigger #;Counts",5,0,5);
151 Add2RawsList(fhRawTrigger ,97, !expert, image, !saveCorr);
153 TH1F* fhRawMean = new TH1F("hRawMean","online mean signal, physics event;", 1000,10000,10000);
154 Add2RawsList( fhRawMean,98, expert, !image, !saveCorr);
155 TH1F* fhRawVertex = new TH1F("hRawVertex","online vertex signal; counts", 100,0,600);
156 Add2RawsList( fhRawVertex,99, expert, !image, !saveCorr);
157 TH1F* fhRawORA = new TH1F("hRawORA","online OR A; counts", 10000,0,10000);
158 Add2RawsList( fhRawORA,100, expert, !image, !saveCorr);
159 TH1F* fhRawORC = new TH1F("hRawORC","online OR C;counts", 10000,0,10000);
160 Add2RawsList( fhRawORC,101, expert, !image, !saveCorr);
162 for (Int_t i=0; i<24; i++)
164 // for events with trigger CALIBRATION_EVENT
165 timeCalname ="hRawCFDcal";
166 ledCalname = "hRawLEDcal";
167 ampCalname = "hRawLEDminCFDcal";
168 qtcCalname = "hRawQTCcal";
173 fhRawCFDcal[i] = new TH1F(timeCalname.Data(), Form("%s;Time [ns];Counts", timeCalname.Data()),10000,0,10000);
174 Add2RawsList( fhRawCFDcal[i],101+i+1, !expert, image, !saveCorr);
175 fhRawLEDcal[i] = new TH1F(ledCalname.Data(), Form("%s;Time [ns];Counts", ledCalname.Data()),10000,0,10000);
176 Add2RawsList( fhRawLEDcal[i],101+i+24+1, !expert, image, !saveCorr);
177 fhRawLEDampcal[i] = new TH1F(ampCalname.Data(), Form("%s;Amplitude [ADC counts];Counts", ampCalname.Data()),1000,0,1000);
178 Add2RawsList( fhRawLEDampcal[i],101+i+48+1, !expert, image, !saveCorr);
179 fhRawQTCcal[i] = new TH1F(qtcCalname.Data(), Form("%s;Charge [??];Counts",qtcCalname.Data()),1000,0,7000);
180 Add2RawsList( fhRawQTCcal[i],101+i+72+1, !expert, image, !saveCorr);
183 TH1F* fhRawTriggerCal = new TH1F("hRawTriggerCal"," laser triggers",6,0,6);
184 Add2RawsList(fhRawTriggerCal ,197 , !expert, image, saveCorr);
186 TH1F* fhRawMeanCal = new TH1F("hRawMeanCal","online mean signal, calibration event",
188 Add2RawsList( fhRawMeanCal,198);
189 TH1F* fhRawVertexCal = new TH1F("hRawVertexCal","online vertex signal, calibration event ",
191 Add2RawsList( fhRawVertexCal,199, expert, !image, !saveCorr);
192 TH1F* fhRawORAcal = new TH1F("hRawORAcal","laser OR A; counts", 10000,0,10000);
193 Add2RawsList( fhRawORAcal,200, expert, !image, !saveCorr );
194 TH1F* fhRawORCcal = new TH1F("hRawORCcal","laserOR C;counts ", 10000,0,10000);
195 Add2RawsList( fhRawORCcal,201, expert, !image, !saveCorr);
196 TH1F* fhMultcal = new TH1F("hMultcal","full mulltiplicity;Multiplicity;Entries", 10000,0,10000);
197 Add2RawsList( fhMultcal,202, expert, !image, !saveCorr );
198 TH1F* fhMultScal = new TH1F("hMultScal","full multiplicity with semi-central trigger;Multiplicity;Entries",
200 Add2RawsList( fhMultScal,203, expert, !image, !saveCorr);
201 TH1F* fhMultCcal = new TH1F("hMultCcal","full multiplicity with central trigger;Multiplicity;Entries",
203 Add2RawsList( fhMultCcal,204, expert, !image, !saveCorr);
205 TH2F* fhEffCFD = new TH2F("hEffCFD","#PMT; #CFD counts/nEvents",24, 0 ,24, 50, 0,5);
206 fhEffCFD->SetOption("COLZ");
207 Add2RawsList( fhEffCFD,205, !expert, !image, saveCorr);
208 TH2F* fhEffLED = new TH2F("hEffLED","#PMT; #LED counts/nEvent",24, 0 ,24,
210 fhEffLED->SetOption("COLZ");
211 Add2RawsList( fhEffLED,206, !expert, !image, saveCorr);
212 TH2F* fhEffQTC = new TH2F("hEffQTC","#PMT; QTC efficiency%s;",24, 0 ,24, 100,0,5);
213 fhEffQTC->SetOption("COLZ");
214 Add2RawsList( fhEffQTC,207, !expert, !image, saveCorr);
216 TH2F* fhCFDcal = new TH2F("hCFDcal","#PMT; CFD {#channnels}",25, 0 ,25, 1000,0,5000);
217 fhCFDcal->SetOption("COLZ");
218 Add2RawsList( fhCFDcal,208, !expert, image, !saveCorr);
219 TH2F* fhLEDcal = new TH2F("hLEDcal","#PMT; LED [#channnels]",25, 0 ,25, 1000,0,5000);
220 fhLEDcal->SetOption("COLZ");
221 Add2RawsList( fhLEDcal,209, !expert, image, !saveCorr);
223 const Char_t *triggers[6] = {"mean", "vertex","ORA","ORC","central","semi-central"};
224 for (Int_t itr=0; itr<6; itr++) {
225 GetRawsData(197)->Fill(triggers[itr], fNumTriggersCal[itr]);
226 GetRawsData(197)->SetBinContent(itr+1, fNumTriggersCal[itr]);
230 //____________________________________________________________________________
231 void AliT0QADataMakerRec::InitDigits()
233 // create Digits histograms in Digits subdir
234 const Bool_t expert = kTRUE ;
235 const Bool_t image = kTRUE ;
237 TH2F * fhDigCFD = new TH2F("fhDigCFD", " CFD digits; #PMT; CFD digits[#channels]",25,-0.5,24.5,100,0,1000);
238 fhDigCFD->SetOption("COLZ");
239 Add2DigitsList( fhDigCFD,0, !expert, image);
240 TH2F *fhDigLEDamp = new TH2F("fhDigLEDamp", " LED-CFD digits; #PMT; LED-CFD amplitude ",25,-0.5,24.5,100,100,1000);
241 fhDigLEDamp->SetOption("COLZ");
242 Add2DigitsList( fhDigLEDamp,1, !expert, !image);
243 TH2F * fhDigQTC = new TH2F("fhDigQTC", " QTC digits; #PMT; QTC amplitude",25,-0.5,24.5,100,100,10000);
244 fhDigQTC->SetOption("COLZ");
245 Add2DigitsList( fhDigQTC,2, !expert, !image);}
247 //____________________________________________________________________________
249 void AliT0QADataMakerRec::InitRecPoints()
251 // create cluster histograms in RecPoint subdir
252 const Bool_t expert = kTRUE ;
253 const Bool_t image = kTRUE ;
255 TH2F* fhRecCFD = new TH2F("hRecCFD"," CFD time;Time [ns];Counts",24, 0 ,24,
257 fhRecCFD->SetOption("COLZ");
258 Add2RecPointsList ( fhRecCFD,0, !expert, image);
260 TH2F* fhRecAmpDiff = new TH2F("hRecAmpDiff"," LED-CFD min QTC amplitude;Amplitude [ADC counts];Counts",
261 24, 0 ,24, 200,-10,10);
262 fhRecAmpDiff->SetOption("COLZ");
263 Add2RecPointsList (fhRecAmpDiff, 1, !expert, image);
265 TH1F *fhMean = new TH1F("hMean","online - rec mean;online - rec mean[#channels];",2000, -1000, 1000);
266 Add2RecPointsList ( fhMean,2, !expert, image);
270 //____________________________________________________________________________
271 void AliT0QADataMakerRec::InitESDs()
273 //create ESDs histograms in ESDs subdir
274 const Bool_t expert = kTRUE ;
275 const Bool_t image = kTRUE ;
277 TH1F *fhESDMean = new TH1F("hESDmean"," ESD mean; mean time[%channels]",1000,0,1000);
278 Add2ESDsList(fhESDMean, 0, !expert, image) ;
279 TH1F * fhESDVertex = new TH1F("hESDvertex","ESDvertex; vertex[cm];",82,-30,30);
280 Add2ESDsList(fhESDVertex, 1, !expert, image) ;
285 //____________________________________________________________________________
286 void AliT0QADataMakerRec::MakeRaws( AliRawReader* rawReader)
290 //fills QA histos for RAW
294 AliT0RawReader *start = new AliT0RawReader(rawReader);
297 AliDebug(AliQAv1::GetQADebugLevel(),Form(" no raw data found!!"));
302 UInt_t type =rawReader->GetType();
303 Int_t allData[110][5];
304 for (Int_t i0=0; i0<105; i0++)
306 for (Int_t j0=0; j0<5; j0++) allData[i0][j0]=0;
308 for (Int_t i=0; i<105; i++)
309 for (Int_t iHit=0; iHit<5; iHit++)
310 allData[i][iHit]= start->GetData(i,iHit);
312 if (allData[0][0]>0) GetRawsData(0) -> Fill( allData[0][0]);
313 Int_t refpoint = allData[0][0];
314 // allData[0][0] = allData[0][0] - 7000;
315 if (type == 8) {shift=101; refpoint=allData[0][0] - 5000;}
316 if (type == 7) shift=0;
318 for (Int_t ik = 0; ik<12; ik++){
319 for (Int_t iHt=0; iHt<1; iHt++){
321 if(allData[ik+1][iHt]>0) {
322 GetRawsData(shift+ik+1) ->
323 Fill(allData[ik+1][iHt]-refpoint);
326 GetRawsData(208)->Fill(ik+1, allData[ik+1][iHt]-refpoint);
329 if(allData[ik+13][iHt] > 0) {
330 GetRawsData(shift+ik+24+1)->
331 Fill(allData[ik+13][iHt]-refpoint);
334 GetRawsData(209)->Fill(ik+1, allData[ik+13][iHt]-refpoint);
339 if(allData[ik+13][iHt] > 0 && allData[ik+1][iHt] >0 )
340 GetRawsData(shift+ik+48+1)->
341 Fill(allData[ik+13][iHt]-allData[ik+1][iHt]);
343 if(allData[2*ik+25][iHt] > 0 || allData[2*ik+26][iHt] > 0) {
344 GetRawsData(shift+ik+72+1)->
345 Fill(allData[2*ik+26][iHt]-allData[2*ik+25][iHt]);
346 if(type == 8 ) feffqtc[ik]++;
350 effic = Float_t(feffC[ik])/Float_t(fnEvent);
351 GetRawsData(205)->Fill(ik+1,effic );
353 effic = Float_t(feffA[ik])/Float_t(fnEvent);
354 GetRawsData(206)->Fill(ik+1,effic );
356 effic = Float_t(feffqtc[ik])/Float_t(fnEvent);
357 GetRawsData(207)->Fill(ik+1, effic);
359 for (Int_t ik = 12; ik<24; ik++) {
360 for (Int_t iHt=0; iHt<1; iHt++) {
361 if(allData[ik+45][iHt]>0) {
363 GetRawsData(shift+ik+1)->
364 Fill(allData[ik+45][iHt]-refpoint);
367 GetRawsData(208)->Fill(ik+1, allData[ik+45][iHt]-refpoint);
371 if(allData[ik+57][iHt] > 0 ) {
372 GetRawsData(shift+ik+24+1)->
373 Fill(allData[ik+57][iHt]-refpoint);
376 GetRawsData(209)->Fill(ik+1, allData[ik+57][iHt]-refpoint);
380 if(allData[2*ik+57][iHt]>0 || allData[2*ik+58][iHt]>0)
382 GetRawsData(shift+ik+72+1)->
383 Fill(allData[2*ik+58][iHt]-allData[2*ik+57][iHt]);
384 if(type == 8 ) feffqtc[ik]++;
387 if(allData[ik+57][iHt] > 0 &&allData[ik+45][iHt]>0)
388 GetRawsData(shift+ik+48+1)->
389 Fill(allData[ik+57][iHt]-allData[ik+45][iHt]);
391 effic = Float_t(feffC[ik])/Float_t(fnEvent);
392 GetRawsData(205)->Fill(ik,effic );
393 effic = Float_t(feffA[ik])/Float_t(fnEvent);
394 GetRawsData(206)->Fill(ik,effic );
396 effic = Float_t(feffqtc[ik])/Float_t(fnEvent);
397 GetRawsData(207)->Fill(ik+1, effic);
400 Int_t trChannel[6] = {49,50,51,52,55,56};
404 for (Int_t iHt=0; iHt<6; iHt++) {
405 for (Int_t itr=0; itr<6; itr++) {
406 if(allData[trChannel[itr]][iHt]>0) fNumTriggers[itr]++;
413 for (Int_t iHt=0; iHt<5; iHt++) {
414 for (Int_t itr=0; itr<6; itr++) {
415 if(allData[trChannel[itr]][iHt]>0)
418 GetRawsData(198+itr)->Fill(allData[trChannel[itr]][iHt]-allData[1][0]);
420 fNumTriggersCal[itr]++;
423 if(allData[53][iHt]>0 && allData[54][iHt]>0)
424 GetRawsData(204)->Fill(allData[53][iHt]-allData[54][iHt]);
434 //____________________________________________________________________________
435 void AliT0QADataMakerRec::MakeDigits( TTree *digitsTree)
437 //fills QA histos for Digits
439 TArrayI *digCFD = new TArrayI(24);
440 TArrayI *digLED = new TArrayI(24);
441 TArrayI *digQT0 = new TArrayI(24);
442 TArrayI *digQT1 = new TArrayI(24);
445 TBranch *brDigits=digitsTree->GetBranch("T0");
446 AliT0digit *fDigits = new AliT0digit() ;
448 brDigits->SetAddress(&fDigits);
450 AliError(Form("EXEC Branch T0 digits not found"));
453 digitsTree->GetEvent(0);
454 digitsTree->GetEntry(0);
455 brDigits->GetEntry(0);
456 fDigits->GetTimeCFD(*digCFD);
457 fDigits->GetTimeLED(*digLED);
458 fDigits->GetQT0(*digQT0);
459 fDigits->GetQT1(*digQT1);
460 refpoint = fDigits->RefPoint();
461 for (Int_t i=0; i<24; i++)
463 if (digCFD->At(i)>0) {
464 Int_t cfd=digCFD->At(i)- refpoint;
465 GetDigitsData(0) ->Fill(i,cfd);
466 GetDigitsData(1) -> Fill(i, (digLED->At(i) - digCFD->At(i)));
467 GetDigitsData(2) -> Fill(i, (digQT1->At(i) - digQT0->At(i)));
478 //____________________________________________________________________________
479 void AliT0QADataMakerRec::MakeRecPoints(TTree * clustersTree)
481 //fills QA histos for clusters
483 AliT0RecPoint* frecpoints= new AliT0RecPoint ();
485 AliError(":MakeRecPoints >> no recpoints found");
488 TBranch *brRec =clustersTree ->GetBranch("T0");
490 brRec->SetAddress(&frecpoints);
492 AliError(Form("EXEC Branch T0 rec not found "));
498 for ( Int_t i=0; i<24; i++) {
500 GetRecPointsData(0) -> Fill(i, frecpoints -> GetTime(i) - frecpoints -> GetTime(0));
502 GetRecPointsData(0) -> Fill(i, frecpoints -> GetTime(i) - frecpoints -> GetTime(12));
503 GetRecPointsData(1) -> Fill( i, frecpoints -> GetAmp(i) - frecpoints->AmpLED(i));
505 Double_t mmm=frecpoints->GetOnlineMean()- frecpoints->GetMeanTime();
506 GetRecPointsData(2) ->Fill(mmm);
510 //____________________________________________________________________________
511 void AliT0QADataMakerRec::MakeESDs(AliESDEvent * esd)
513 //fills QA histos for ESD
515 GetESDsData(0) -> Fill(esd->GetT0());
516 GetESDsData(1)-> Fill(esd->GetT0zVertex());