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 ---
35 #include "AliESDEvent.h"
37 #include "AliT0digit.h"
39 #include "AliT0RecPoint.h"
40 #include "AliT0QADataMakerRec.h"
41 #include "AliQAChecker.h"
42 #include "AliT0RawReader.h"
43 #include "AliT0RecoParam.h"
44 #include "THnSparse.h"
46 #include "Riostream.h"
47 ClassImp(AliT0QADataMakerRec)
49 //____________________________________________________________________________
50 AliT0QADataMakerRec::AliT0QADataMakerRec() :
51 AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kT0),
52 "T0 Quality Assurance Data Maker")
57 // RS: There is some inconsistency here: the separation of physics and calib. events/histos is done by
58 // fEventSpecie. Why do we book separate histos on different slots for calib and physics ?
59 // I am changing this in such way that we don't need local counters like fNumTriggers (the corresponding
60 // histos now incremented in the MakeRaws, and for the normalization I will use the framework's counters
61 // AliQADataMaker::GetEvCountCycle(...), AliQADataMaker::GetEvCountTotal(...)
62 // All these fTrEff.. feff.. will by directly filled in corresponding histos
64 for(Int_t i=0; i<24; i++) fMeans[i]=2500;
68 //____________________________________________________________________________
69 AliT0QADataMakerRec::AliT0QADataMakerRec(const AliT0QADataMakerRec& qadm) :
74 SetName((const char*)qadm.GetName()) ;
75 SetTitle((const char*)qadm.GetTitle());
78 //__________________________________________________________________
79 AliT0QADataMakerRec& AliT0QADataMakerRec::operator = (const AliT0QADataMakerRec& qadm )
82 this->~AliT0QADataMakerRec();
83 new(this) AliT0QADataMakerRec(qadm);
86 //__________________________________________________________________
87 AliT0QADataMakerRec::~AliT0QADataMakerRec()
91 //____________________________________________________________________________
92 void AliT0QADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
94 //Detector specific actions at end of cycle
96 AliInfo(Form("Task: %d",task));
97 ResetEventTrigClasses();
98 AliQAChecker::Instance()->Run(AliQAv1::kT0, task, list) ;
99 // const Char_t *triggers[6] = {"mean", "vertex","ORA","ORC","central","semi-central"};
102 for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
104 // RS: There is some inconsistency here: the separation of physics and calib. events/histos is done by
105 // fEventSpecie. Why do we book separate histos on different slots for calib and physics ?
106 // I am changing this in such way that we don't need local counters like fNumTriggers (the corresponding
107 // histos now incremented in the MakeRaws, and for the normalization I will use the framework's counters
108 // AliQADataMaker::GetEvCountCycle(...), AliQADataMaker::GetEvCountTotal(...)
110 // I think the histos xx+250 should be suppressed (the xx calib histos of specie==calibration will be
111 // used automatically)
113 if (! IsValidEventSpecie(specie, list)) continue;
114 SetEventSpecie(AliRecoParam::ConvertIndex(specie));
116 for (int itc=-1;itc<GetNTrigClasses();itc++) { // RS: loop over eventual clones per trigger class
118 if ( task == AliQAv1::kRAWS ) {
120 float nEvent = GetEvCountCycleRaws(itc); // counted events for given trigger class
122 if ( (htmp=GetRawsData(169,itc)) ) htmp->Scale(1./nEvent); // RS
123 if ( (htmp=GetRawsData(207,itc)) ) htmp->Scale(1./nEvent); // RS // former feffPhysC count.
124 if ( (htmp=GetRawsData(208,itc)) ) htmp->Scale(1./nEvent); // RS // former feffA count. Remove?
125 if ( (htmp=GetRawsData(209,itc)) ) htmp->Scale(1./nEvent); // RS // former feffqtc count.
128 } // RS: loop over eventual clones per trigger class
129 } // loop over species
134 //____________________________________________________________________________
135 void AliT0QADataMakerRec::StartOfDetectorCycle()
137 //Detector specific actions at start of cycle
141 //____________________________________________________________________________
142 void AliT0QADataMakerRec::InitRaws()
144 // create Raw histograms in Raw subdir
145 const Bool_t expert = kTRUE ;
146 const Bool_t saveCorr = kTRUE ;
147 const Bool_t image = kTRUE ;
151 const Char_t *triggers[6] = {"mean", "vertex","ORA","ORC","central","semi-central"};
154 for (Int_t i=0; i<500; i++){
160 TString timename, ampname, qtcname, ledname;
161 TString timeCalname, ampCalname, ledCalname, qtcCalname;
162 TString qt1name, qt0name, qt1Calname, qt0Calname;
165 TH1F* fhRefPoint = new TH1F("hRefPoint","Ref Point", 10000, 0 ,50000);
166 fhRefPoint->SetLabelSize(0.02);
167 Add2RawsList( fhRefPoint,0, expert, !image, !saveCorr);
171 TH1F * fhRawLEDamp[24];
172 TH1F *fhRawQTC[24]; TH1F * fhRawLED[24];
173 TH1F *fhRawQT1[24]; TH1F *fhRawQT0[24];
174 TH1F* fhRawNhits[24];
176 for (Int_t i=0; i<24; i++)
181 qt0name = "hRawQT0_";
182 qt1name = "hRawQT1_";
183 ampname = "hRawLEDminCFD";
192 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]);
193 // ForbidCloning(fhRawCFD[i]); //RS I don't know how histos 1-24 should be processed in MakeRaws, for the moment forbidding the cloning
194 Add2RawsList( fhRawCFD[i],i+1, expert, !image, !saveCorr);
195 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]);
196 Add2RawsList( fhRawLED[i],i+25, expert, !image, !saveCorr);
197 fhRawLEDamp[i] = new TH1F(ampname.Data(), Form("%s;LED-CFD [#channels];Counts", ampname.Data()),1000,0,1000);
198 Add2RawsList( fhRawLEDamp[i],i+49, expert, !image, !saveCorr);
199 fhRawQTC[i] = new TH1F(qtcname.Data(), Form("%s;QTC[#channels];Counts", qtcname.Data()),10000,0,10000);
200 Add2RawsList( fhRawQTC[i],i+73, expert, !image, !saveCorr);
201 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]);
202 Add2RawsList( fhRawQT1[i],97+i, expert, !image, !saveCorr);
203 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]);
204 Add2RawsList( fhRawQT0[i],121+i, expert, !image, !saveCorr);
206 fhRawNhits[i] = new TH1F(nhits.Data(), Form("%s;#Hits;Events", nhits.Data()),20, 0, 20);
207 Add2RawsList( fhRawNhits[i],176+i, expert, !image, !saveCorr);
212 TH1F* fhRawTrigger = new TH1F("hRawTrigger"," triggers;Trigger ;Counts",6,0,6);
213 for (Int_t itr=0; itr<6; itr++) fhRawTrigger->Fill(triggers[itr], 0); // RS Modified to allow cloning (no fNumTriggers member anymore)
214 fhRawTrigger->SetMinimum(0);
215 fhRawTrigger->SetMaximum(2);
216 Add2RawsList(fhRawTrigger ,169, !expert, image, !saveCorr);
218 TH1F* fhRawMean = new TH1F("hRawMean","online mean signal, physics event;",Int_t((high[170]-low[170])/4),low[170],high[170]);
219 Add2RawsList( fhRawMean,170, expert, !image, !saveCorr);
221 TH1F* fhRawVertex = new TH1F("hRawVertex","online vertex signal; counts",Int_t((high[171]-low[171])/4),low[171],high[171]);
222 Add2RawsList( fhRawVertex,171, expert, image, !saveCorr);
224 TH1F* fhRawORA = new TH1F("hRawORA","online OR A; counts",Int_t((high[172]-low[172])/4),low[172],high[172]);
225 Add2RawsList( fhRawORA,172, expert, !image, !saveCorr);
226 TH1F* fhRawORC = new TH1F("hRawORC","online OR C;counts",Int_t(( high[173]-low[173])/4),low[173],high[173]);
227 Add2RawsList( fhRawORC,173, expert, !image, !saveCorr);
228 TH1F* fhMultCentr = new TH1F("hMultCentr","online trigger Central;counts ",Int_t(( high[174]-low[174])/4),low[174],high[174]);
229 Add2RawsList( fhMultCentr,174, expert, !image, !saveCorr);
230 TH1F* fhMultSeCentr = new TH1F("hMultSemiCentr","online trigger SemiCentral;counts ",Int_t(( high[175]-low[175])/4),low[175],high[175]);
231 Add2RawsList( fhMultSeCentr,175, expert, !image, !saveCorr);
233 TH1F* fhMultA = new TH1F("hMultA","full mulltiplicity A side;Multiplicity;Entries", Int_t((high[201]-low[201])/4) ,low[201],high[201]);
234 Add2RawsList( fhMultA,201, expert, image, !saveCorr );
236 TH1F* fhMultAS = new TH1F("hMultASemi","full multiplicity with semi-central trigger A side ;Multiplicity;Entries",
237 Int_t((high[202]-low[202])/4),low[202],high[202] );
238 Add2RawsList( fhMultAS, 202, expert, !image, !saveCorr);
239 TH1F* fhMultAC = new TH1F("hMultACentr","full multiplicity with central trigger;Multiplicity;Entries",
240 Int_t((high[203]-low[203])/4),low[203],high[203]);
241 Add2RawsList( fhMultAC, 203, expert, !image, !saveCorr);
245 TH1F* fhMultC = new TH1F("hMultC","full mulltiplicity C side;Multiplicity;Entries", Int_t(high[204]-low[204]/4) ,low[204],high[204]);
246 Add2RawsList( fhMultC,204, expert, image, !saveCorr );
247 TH1F* fhMultCS = new TH1F("hMultCSemi","full multiplicity with semi-central trigger C side;Multiplicity;Entries",
248 Int_t((high[205]-low[205])/4),low[205],high[205] );
249 Add2RawsList( fhMultCS,205, expert, !image, !saveCorr);
250 TH1F* fhMultCC = new TH1F("hMultCCentr","full multiplicity with central trigger C side;Multiplicity;Entries",
251 Int_t((high[206]-low[206])/4),low[206],high[206]);
252 Add2RawsList( fhMultCC,206, expert, !image, !saveCorr);
256 TH1F* fhCFDeff= new TH1F("fhCFDeff"," CFD efficiency; #PMT; #CFD counts/nEvents",24, 0 ,24);
257 fhCFDeff->SetMinimum(0);
258 fhCFDeff->SetMaximum(2);
259 Add2RawsList( fhCFDeff, 207, !expert, image, !saveCorr);
260 TH1F* fhEffLED = new TH1F("hEffLED","LEDefficiecy; #PMT; #LED counts/nEvent",24, 0 ,24);
261 fhEffLED ->SetMinimum(0);
262 fhEffLED->SetMaximum(2);
263 Add2RawsList( fhEffLED,208, !expert, image, !saveCorr);
265 TH1F* fhEffQTC = new TH1F("hEffQTC","QTC efficiecy; #PMT; QTC efficiency%s;",24, 0 ,24);
266 fhEffQTC->SetMinimum(0);
267 fhEffQTC->SetMaximum(2);
268 Add2RawsList( fhEffQTC,209, !expert, image, !saveCorr);
270 TH2F* fhCFD = new TH2F("hCFD","CFD ; #PMT; CFD {#channnels}",25, 0 ,25,Int_t((high[210]-low[210])/4),low[210],high[210]);
271 fhCFD->SetOption("COLZ");
272 Add2RawsList( fhCFD,210, expert, image, !saveCorr);
274 TH2F* fhLED = new TH2F("hLED","LED ; #PMT; LED [#channnels]",25, 0 ,25,Int_t((high[211]-low[211])/4),low[211],high[211]);
275 fhLED->SetOption("COLZ");
276 Add2RawsList( fhLED,211, expert, image, !saveCorr);
278 TH2F* fhQTC = new TH2F("hQTC","QTC ; #PMT; QTC [#channnels]",25, 0 ,25,Int_t( high[212]-low[212]),low[212],high[212]);
279 fhQTC->SetOption("COLZ");
280 Add2RawsList( fhQTC,212, expert, image, !saveCorr);
282 TH1F* fhNumPMTA= new TH1F("hNumPMTA","number of PMT hitted per event A side",13, 0 ,13);
283 Add2RawsList(fhNumPMTA ,213, expert, image, !saveCorr);
285 TH1F* fhNumPMTC= new TH1F("hNumPMTC","number of PMT hitted per event C side",13, 0 ,13);
286 Add2RawsList(fhNumPMTC ,214, expert, image, !saveCorr);
288 TH1F* fhHitsOrA= new TH1F("fhHitsOrA","T0_OR A hit multiplicitie",20, 0 ,20);
289 Add2RawsList( fhHitsOrA,215, expert, !image, !saveCorr);
291 TH1F* fhHitsOrC= new TH1F("fhHitsOrC","T0_OR C hit multiplicitie",20, 0 ,20);
292 Add2RawsList(fhHitsOrC ,216, expert, !image, !saveCorr);
295 TH1F* fhOrCminOrA= new TH1F("fhOrCminOrA","T0_OR C - T0_OR A",10000,-5000,5000);
296 Add2RawsList( fhOrCminOrA,219, expert, !image, !saveCorr);
298 TH1F* fhOrCminOrATvdcOn= new TH1F("fhOrCminOrATvdcOn","T0_OR C - T0_OR A TVDC on",10000,-5000,5000);
299 Add2RawsList( fhOrCminOrATvdcOn,217, expert, !image, !saveCorr);
302 TH1F* fhOrCminOrATvdcOff= new TH1F("fhOrCminOrATvdcOff","T0_OR C - T0_OR A TVDC off",10000,-5000,5000);
303 Add2RawsList( fhOrCminOrATvdcOff,218, expert, !image, !saveCorr);
305 //satellite & beam background
306 TH2F* fhBeam = new TH2F("fhBeam", " Mean vs Vertex ", 120, -30, 30, 120, -30, 30);
307 fhBeam->SetOption("COLZ");
308 fhBeam->GetXaxis()->SetTitle("(T0C-T0A/2, ns");
309 fhBeam->GetYaxis()->SetTitle("(T0C+T0A/2, ns");
310 Add2RawsList( fhBeam,220, !expert, image, !saveCorr);
311 TH2F* fhBeamTVDCon = new TH2F("fhBeamTVDCon", " Mean vs Vertex TVDC on ",50, -5, 5, 50, -5, 5) ;
312 fhBeamTVDCon->SetOption("COLZ");
313 fhBeamTVDCon->GetXaxis()->SetTitle("(T0C-T0A/2, ns");
314 fhBeamTVDCon->GetYaxis()->SetTitle("(T0C+T0A/2, ns");
315 Add2RawsList( fhBeamTVDCon,221, expert, image, !saveCorr);
316 TH2F* fhBeamTVDCoff = new TH2F("fhBeamTVDCoff", " Mean vs Vertex TVDC off", 120, -30, 30, 120, -30, 30);
317 fhBeamTVDCoff->GetXaxis()->SetTitle("(T0C-T0A/2, ns");
318 fhBeamTVDCoff->GetYaxis()->SetTitle("(T0C+T0A/2, ns");
319 fhBeamTVDCoff->SetOption("COLZ");
320 Add2RawsList( fhBeamTVDCoff,222, expert, image, !saveCorr);
321 TH1F* fhMean = new TH1F("fhMean", " (T0A+T0C)/2, ps ", 200, -2000, 2000);
322 Add2RawsList( fhMean,223, !expert, image, !saveCorr);
324 TH1F* fhVertex1st = new TH1F("fhVertex1st", " (T0A-T0C)/2, ps 1st particle", 200, -2000, 2000);
325 Add2RawsList( fhVertex1st ,225, !expert, image, !saveCorr);
326 TH1F* fhVertexBest = new TH1F("fhVertexBest", " (T0A-T0C)/2, ps , best particle", 200, -2000, 2000);
327 Add2RawsList( fhVertexBest ,226, !expert, image, !saveCorr);
328 TH1F* fhMean1st = new TH1F("fhMean1st", " (T0A + T0C)/2, ps , 1st particle", 200, -2000, 2000);
329 Add2RawsList( fhMean1st ,227, expert, image, !saveCorr);
332 TH2F* fhBCID = new TH2F("fhBCID", "header BCID vs TRM BC ID ", 500, 0, 5000, 500, 0, 5000);
333 fhBCID->SetOption("COLZ");
334 fhBCID->GetXaxis()->SetTitle("TRM BC ID");
335 fhBCID->GetYaxis()->SetTitle("event header BC ID");
336 Add2RawsList(fhBCID ,224, !expert, image, !saveCorr);
338 ClonePerTrigClass(AliQAv1::kRAWS); // this should be the last line
341 //____________________________________________________________________________
342 void AliT0QADataMakerRec::InitDigits()
344 // create Digits histograms in Digits subdir
345 const Bool_t expert = kTRUE ;
346 const Bool_t image = kTRUE ;
348 TH2F * fhDigCFD = new TH2F("fhDigCFD", " CFD digits; #PMT; CFD digits[#channels]",25,-0.5,24.5,100,0,1000);
349 fhDigCFD->SetOption("COLZ");
350 Add2DigitsList( fhDigCFD,0, !expert, image);
351 TH2F *fhDigLEDamp = new TH2F("fhDigLEDamp", " LED-CFD digits; #PMT; LED-CFD amplitude ",25,-0.5,24.5,100,100,1000);
352 fhDigLEDamp->SetOption("COLZ");
353 Add2DigitsList( fhDigLEDamp,1, !expert, !image);
354 TH2F * fhDigQTC = new TH2F("fhDigQTC", " QTC digits; #PMT; QTC amplitude",25,-0.5,24.5,100,100,10000);
355 fhDigQTC->SetOption("COLZ");
356 Add2DigitsList( fhDigQTC,2, !expert, !image);
358 ClonePerTrigClass(AliQAv1::kDIGITS); // this should be the last line
361 //____________________________________________________________________________
363 void AliT0QADataMakerRec::InitRecPoints()
365 // create cluster histograms in RecPoint subdir
366 const Bool_t expert = kTRUE ;
367 const Bool_t image = kTRUE ;
369 TH2F* fhRecCFD = new TH2F("hRecCFD"," CFD time;#PMT; CFD Time [ns];",24, 0 ,24,
371 fhRecCFD->SetOption("COLZ");
372 Add2RecPointsList ( fhRecCFD,0, !expert, image);
374 TH2F* fhRecAmpDiff = new TH2F("hRecAmpDiff"," LED-CFD min QTC amplitude;#PMT; difference [MIPs];",
375 24, 0 ,24, 200,-10,10);
376 fhRecAmpDiff->SetOption("COLZ");
377 Add2RecPointsList (fhRecAmpDiff, 1, !expert, image);
379 TH1F *fhMean = new TH1F("hMean","online - rec mean;online - rec mean[#channels];",2000, -1000, 1000);
380 Add2RecPointsList ( fhMean,2, !expert, image);
382 ClonePerTrigClass(AliQAv1::kRECPOINTS); // this should be the last line
385 //____________________________________________________________________________
386 void AliT0QADataMakerRec::InitESDs()
388 //create ESDs histograms in ESDs subdir
389 const Bool_t expert = kTRUE ;
390 const Bool_t image = kTRUE ;
392 TH1F *fhESDMean = new TH1F("hESDmean"," ESD mean; mean time[%channels]",1000, -5, 5);
393 Add2ESDsList(fhESDMean, 0, expert, !image) ;
394 TH1F * fhESDVertex = new TH1F("hESDvertex","ESDvertex; vertex[cm];",82,-30,30);
395 Add2ESDsList(fhESDVertex, 1, expert, !image) ;
397 TH1F * fhESDResolution = new TH1F("hESDResolution","(T0A-T0C)/2 corrected by SPD vertex; ns",800,-2,2);
398 Add2ESDsList(fhESDResolution, 2, !expert, image) ;
400 ClonePerTrigClass(AliQAv1::kESDS); // this should be the last line
403 //____________________________________________________________________________
404 void AliT0QADataMakerRec::MakeRaws( AliRawReader* rawReader)
407 for(Int_t i=0; i<24; i++) time[i] = 0;
410 //fills QA histos for RAW
412 // Int_t refPointParam = GetRecoParam()->GetRefPoint();
414 Int_t refPointParam = 0;
416 AliT0RawReader *start = new AliT0RawReader(rawReader);
418 if (! start->Next()) {
419 AliDebug(AliQAv1::GetQADebugLevel(),Form(" no raw data found!!"));
423 UInt_t type =rawReader->GetType();
424 if (GetEventSpecie()==AliRecoParam::kCalib && type!=8) {
429 // RS: Don't use custom counters, they create problems with trigger cloning
430 // Use instead framework counters, incremented in the end of this routine
431 // RS: There is some inconsistency here: the separation of physics and calib. events/histos is done by
432 // fEventSpecie. Why do we book separate histos on different slots for calib and physics ?
433 // I am changing this in such way that we don't need local counters like fNumTriggers (the corresponding
434 // histos now incremented in the MakeRaws, and for the normalization I will use the framework's counters
435 // AliQADataMaker::GetEvCountCycle(...), AliQADataMaker::GetEvCountTotal(...)
437 // I think the histos xx+250 should be suppressed (the xx calib histos of specie==calibration will be
438 // used automatically)
443 UInt_t bcid = rawReader->GetBCID();
444 UInt_t trmbcid = start->GetTRMBunchID();
445 FillRawsData(224,trmbcid, bcid);
447 // if (type == 7){ shift=1; fnEventPhys++;}
448 Int_t allData[110][5];
449 for (Int_t i0=0; i0<110; i0++)
451 for (Int_t j0=0; j0<5; j0++) allData[i0][j0]=0;
453 for (Int_t i=0; i<107; i++)
454 for (Int_t iHit=0; iHit<5; iHit++)
455 allData[i][iHit]= start->GetData(i,iHit);
457 //if ( allData[0][0] > 0 && (type == 7))
458 FillRawsData(0, allData[0][0]);
459 refpoint = allData[refPointParam][0];
460 if (refPointParam < 0 ) refpoint=0;
461 if (refPointParam == 0 ) refpoint = allData[0][0] - 5000;
463 Int_t sideshift, sideshiftqtc;
466 for (Int_t ik = 0; ik<24; ik++) {
470 if(allData[ik+sideshift][0]>0 /*&& type == 7 */ ) numPmtC++;
473 if(allData[ik+45][0]>0 /*&& type == 7 */) numPmtA++;
479 for (Int_t iHt=0; iHt<5; iHt++) {
481 if(allData[ik+sideshift][iHt]>0) {
482 FillRawsData(shift+ik+1, allData[ik+sideshift][iHt]);
483 FillRawsData(210+shift, ik+1, allData[ik+sideshift][iHt]);
484 FillRawsData(207+shift,ik,1.); // instead of incrementing former feff's, increment histo directly
485 AliDebug(50,Form("%i CFD %i data %s",ik, ik+sideshift, GetRawsData(shift+ik+1)->GetName()));
490 if(allData[ik+12+sideshift][iHt] > 0) {
491 FillRawsData(shift+ik+24+1,allData[ik+12+sideshift][iHt]);
492 FillRawsData(211+shift,ik+1, allData[ik+12+sideshift][iHt]);
493 AliDebug(50,Form("%i LED %i data %s",ik, ik+12+sideshift, GetRawsData(shift+ik+1+24)->GetName()));
494 FillRawsData(208+shift,ik,1.); // instead of incrementing former feff's, increment histo directly
499 if(allData[ik+12+sideshift][iHt] > 0 && allData[ik+sideshift][iHt] >0 )
500 FillRawsData(shift+ik+48+1, allData[ik+12+sideshift][iHt]-allData[ik+sideshift][iHt]);
503 if(allData[2*ik+sideshiftqtc+24][iHt] > 0 &&
504 allData[2*ik+sideshiftqtc+25][iHt] > 0) {
505 FillRawsData(shift+ik+72+1, allData[2*ik+sideshiftqtc+24][iHt]-allData[2*ik+sideshiftqtc+25][iHt]);
506 FillRawsData(212+shift,ik+1, allData[2*ik+sideshiftqtc+24][iHt]-allData[2*ik+sideshiftqtc+25][iHt]);
508 FillRawsData(209+shift,ik,1.); // instead of incrementing former feff's, increment histo directly
510 AliDebug(50,Form("%i QTC %i data %s",ik, 2*ik+sideshiftqtc+24, GetRawsData(shift+ik+1+72)->GetName()));
513 if(allData[2*ik+sideshiftqtc+24][iHt] > 0) {
514 AliDebug(50,Form("%i QT0 %i data %s",ik, 2*ik+sideshiftqtc+24, GetRawsData(shift+ik+1+96)->GetName()));
515 FillRawsData(shift+ik+96+1,allData[2*ik+sideshiftqtc+24][iHt]);
517 if(allData[2*ik+sideshiftqtc+25][iHt] > 0) {
518 AliDebug(50,Form("%i QT0 %i data %s",ik, 2*ik+sideshiftqtc+25, GetRawsData(shift+ik+1+120)->GetName()));
519 FillRawsData(shift+ik+120+1,allData[2*ik+sideshiftqtc+25][iHt]);
523 FillRawsData(ik+176,nhitsPMT);
524 FillRawsData(213,numPmtA);
525 FillRawsData(214,numPmtC);
528 Int_t trChannel[6] = {49,50,51,52,55,56};
529 Float_t ch2cm = 24.4*0.029979;
532 for (Int_t iHt=0; iHt<5; iHt++) {
534 //orA-orC phys tvdc 1
535 if((allData[51][iHt]>0 && allData[52][iHt]>0) && allData[50][iHt]>0) {
536 AliDebug(10,Form("orA-orC phys tvdc 1 %i data %s", 217+shift, GetRawsData(shift+217)->GetName()));
538 FillRawsData(217+shift,(allData[52][iHt] - allData[51][iHt])*ch2cm);
540 //orA-orC phys tvdc 0
541 if((allData[51][iHt]>0 && allData[52][iHt]>0) && allData[50][iHt]<=0) {
542 AliDebug(10,Form("orA-orC phys tvdc 0 %i data %s", 218+shift, GetRawsData(shift+218)->GetName()));
544 FillRawsData(218+shift,(allData[52][iHt] - allData[51][iHt])*ch2cm);
546 if(allData[51][iHt]>0 && allData[52][iHt]>0) {
547 AliDebug(50,Form("orA-orC phys tvdc all %i data %s", 219+shift, GetRawsData(shift+219)->GetName()));
548 FillRawsData(219+shift,(allData[52][iHt] - allData[51][iHt])*ch2cm);
550 for (Int_t itr=0; itr<6; itr++) {
551 if (allData[trChannel[itr]][iHt] >0) {
553 // RS instead of incremented custom counters, fill directly the specie-specific histos
554 // FillRawsData(169+shift, 0.5+itr, 1.); // RS: increment counters
555 FillRawsData(169+shift, itr, 1.); // RS: increment counters
556 AliDebug(50,Form(" triggers %i data %s", 170+itr+shift, GetRawsData(170+itr+shift)->GetName()));
558 FillRawsData(170+itr+shift,allData[trChannel[itr]][iHt]);
562 /* if(type == 7) */if(allData[51][iHt] >0) nhitsOrA++;
563 /* if(type == 7) */if(allData[52][iHt] >0) nhitsOrC++;
565 //mult trigger signals phys
567 if(allData[53][iHt]>0 && allData[54][iHt]>0) {
568 AliDebug(50,Form(" mpdA %i data %s", 201+shift, GetRawsData(201+shift)->GetName()));
570 FillRawsData(201+shift,allData[53][iHt]-allData[54][iHt]);
571 if(allData[56][iHt]>0) FillRawsData(202+shift,allData[53][iHt]-allData[54][iHt]);
572 if(allData[55][iHt]>0) FillRawsData(203+shift,allData[53][iHt]-allData[54][iHt]);
576 if(allData[105][iHt]>0 && allData[106][iHt]>0) {
577 AliDebug(50,Form(" mpdC %i data %s", 204+shift, GetRawsData(204+shift)->GetName()));
579 FillRawsData(204+shift,allData[105][iHt]-allData[106][iHt]);
580 if(allData[56][iHt]>0) FillRawsData(205+shift,allData[105][iHt]-allData[106][iHt]);
581 if(allData[55][iHt]>0) FillRawsData(206+shift,allData[105][iHt]-allData[106][iHt]);
585 FillRawsData(215,nhitsOrA);
586 FillRawsData(216,nhitsOrC);
589 for (int itr=-1;itr<GetNEventTrigClasses();itr++) { //RS loop over all active trigger classes, including the global one
590 int itrID = itr==-1 ? -1 : GetEventTrigClass(itr)->GetUniqueID();
591 int nEvent = GetEvCountCycleRaws(itrID);
593 // if(type == 7) { // RS Do we need here event type check, specie should do the job
595 for (Int_t ik=0; ik<24; ik++) {
597 TH1* hik = (TH1*) GetRawsData(ik+1,itrID); if (!hik) continue;
598 if ( hik -> GetEntries() < 100) continue;
599 hik->GetXaxis()->SetRangeUser(2000, 3000);
600 int maxBin = hik->GetMaximumBin();
601 double meanEstimate = hik->GetBinCenter( maxBin);
602 TF1* fit= new TF1("fit","gaus", meanEstimate - 40, meanEstimate + 40);
603 fit->SetParameters (hik->GetBinContent(maxBin), meanEstimate, 80);
604 hik->Fit("fit","RQ","Q", meanEstimate-40, meanEstimate+40);
605 fMeans[ik]= (Int_t) fit->GetParameter(1);
606 hik->GetXaxis()->SetRangeUser(0, 30000);
611 TH2 *h220 = (TH2*)GetRawsData(220,itrID);
612 TH2 *h221 = (TH2*)GetRawsData(221,itrID);
613 TH2 *h222 = (TH2*)GetRawsData(222,itrID);
614 TH1 *h223 = (TH1*)GetRawsData(223,itrID);
615 TH1 *h225 = (TH1*)GetRawsData(225,itrID);
616 TH1 *h226 = (TH1*)GetRawsData(226,itrID);
617 TH1 *h227 = (TH1*)GetRawsData(227,itrID);
618 if( nEvent>1000 && (h220 || h221 || h222 || h223 || h225 || h226|| h227) )
620 Int_t besttimeA=9999999;
621 Int_t besttimeC=9999999;
622 Int_t time1stA=9999999;
623 Int_t time1stC=9999999;
624 for (Int_t ipmt=0; ipmt<12; ipmt++){
625 if(allData[ipmt+1][0] > 1 ) {
626 // time[ipmt] = allData[ipmt+1][0] - Int_t(GetRawsData(1)->GetMean());
627 time[ipmt] = allData[ipmt+1][0] - fMeans[ipmt];
628 if(TMath::Abs(time[ipmt]) < TMath::Abs(besttimeC))
629 besttimeC=time[ipmt]; //timeC
630 if(time[ipmt] < time1stC) time1stC=time[ipmt]; //timeC
633 for ( Int_t ipmt=12; ipmt<24; ipmt++){
634 if(allData[ipmt+45][0] > 0) {
635 time[ipmt] = allData[ipmt+45][0] - fMeans[ipmt] ;
636 if(TMath::Abs(time[ipmt]) < TMath::Abs(besttimeA) )
637 besttimeA=time[ipmt]; //timeA
638 if(time[ipmt] < time1stA) time1stA=time[ipmt]; //timeC
641 if(besttimeA<99999 &&besttimeC< 99999) {
642 Float_t t0 = 24.4 * (Float_t( besttimeA+besttimeC)/2. );
643 Float_t ver = 24.4 * Float_t( besttimeC-besttimeA)/2.;
644 if (h220) h220->Fill(0.001*ver, 0.001*(t0));
645 if(allData[50][0] > 0) if (h221) h221->Fill(0.001*ver, 0.001*(t0));
646 if(allData[50][0] <= 0) if (h222) h222->Fill(0.001*ver, 0.001*(t0));
647 if (h223) h223->Fill(t0);
648 if (h226) h226->Fill(ver);
650 if(time1stA<99999 &&time1stC< 99999) {
651 Float_t t01st = 24.4 * (Float_t( time1stA + time1stC)/2. );
652 Float_t ver1st = 24.4 * Float_t( time1stC - time1stA)/2.;
653 if (h225) h225->Fill(ver1st);
654 if (h227) h227->Fill(t01st);
658 } // RS loop over all active trigger classes, including the global one
660 IncEvCountCycleRaws();
661 IncEvCountTotalRaws();
666 //____________________________________________________________________________
667 void AliT0QADataMakerRec::MakeDigits( TTree *digitsTree)
669 //fills QA histos for Digits
671 TArrayI *digCFD = new TArrayI(24);
672 TArrayI *digLED = new TArrayI(24);
673 TArrayI *digQT0 = new TArrayI(24);
674 TArrayI *digQT1 = new TArrayI(24);
677 TBranch *brDigits=digitsTree->GetBranch("T0");
678 AliT0digit *fDigits = new AliT0digit() ;
680 brDigits->SetAddress(&fDigits);
682 AliError(Form("EXEC Branch T0 digits not found"));
689 digitsTree->GetEvent(0);
690 digitsTree->GetEntry(0);
691 brDigits->GetEntry(0);
692 fDigits->GetTimeCFD(*digCFD);
693 fDigits->GetTimeLED(*digLED);
694 fDigits->GetQT0(*digQT0);
695 fDigits->GetQT1(*digQT1);
696 refpoint = fDigits->RefPoint();
697 for (Int_t i=0; i<24; i++) {
698 if (digCFD->At(i)>0) {
699 Int_t cfd=digCFD->At(i)- refpoint;
700 FillDigitsData(0,i,cfd);
701 FillDigitsData(1,i, (digLED->At(i) - digCFD->At(i)));
702 FillDigitsData(2,i, (digQT1->At(i) - digQT0->At(i)));
711 IncEvCountCycleDigits();
712 IncEvCountTotalDigits();
717 //____________________________________________________________________________
718 void AliT0QADataMakerRec::MakeRecPoints(TTree * clustersTree)
720 //fills QA histos for clusters
722 AliT0RecPoint* frecpoints= new AliT0RecPoint ();
724 AliError(":MakeRecPoints >> no recpoints found");
727 TBranch *brRec =clustersTree ->GetBranch("T0");
729 brRec->SetAddress(&frecpoints);
731 AliError(Form("EXEC Branch T0 rec not found "));
737 for ( Int_t i=0; i<24; i++) {
739 FillRecPointsData(0, i, frecpoints -> GetTime(i) - frecpoints -> GetTime(0));
741 FillRecPointsData(0, i, frecpoints -> GetTime(i) - frecpoints -> GetTime(12));
742 FillRecPointsData(1, i, frecpoints -> GetAmp(i) - frecpoints->AmpLED(i));
744 Double_t mmm=frecpoints->GetOnlineMean()- frecpoints->GetMeanTime();
745 FillRecPointsData(2,mmm);
747 IncEvCountCycleRecPoints();
748 IncEvCountTotalRecPoints();
753 //____________________________________________________________________________
754 void AliT0QADataMakerRec::MakeESDs(AliESDEvent * esd)
756 //fills QA histos for ESD
758 const Double32_t *mean;
759 mean = esd->GetT0TOF();
760 Double32_t t0time= 0.001*mean[0];
761 Double32_t orA= 0.001*mean[1];
762 Double32_t orC=0.001* mean[2];
764 if (t0time<99) FillESDsData(0,t0time);
765 if( esd->GetT0zVertex() <99) FillESDsData(1, esd->GetT0zVertex());
766 if( orA<99 && orC<99) FillESDsData(2,(orA-orC)/2.);
768 IncEvCountCycleESDs();
769 IncEvCountTotalESDs();
772 //____________________________________________________________________________
773 //____________________________________________________________________________
774 void AliT0QADataMakerRec::ResetDetector(AliQAv1::TASKINDEX_t task)
777 //reset the detector histograms for a given task
778 AliQADataMakerRec::ResetDetector(task);
783 void AliT0QADataMakerRec::GetMeanAndSigma(TH1F* hist, Float_t &mean, Float_t &sigma)
786 const double window = 3.; //fit window
788 double meanEstimate, sigmaEstimate;
790 maxBin = hist->GetMaximumBin(); //position of maximum
791 meanEstimate = hist->GetBinCenter( maxBin); // mean of gaussian sitting in maximum
792 sigmaEstimate = hist->GetRMS();
793 TF1* fit= new TF1("fit","gaus", meanEstimate - window*sigmaEstimate, meanEstimate + window*sigmaEstimate);
794 fit->SetParameters(hist->GetBinContent(maxBin), meanEstimate, sigmaEstimate);
795 hist->Fit("fit","RQ","Q");
797 mean = (Float_t) fit->GetParameter(1);
798 sigma = (Float_t) fit->GetParameter(2);