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 **************************************************************************/
16 ///////////////////////////////////////////////////////////////////////
18 // Produces the data needed to calculate the TOF quality assurance. //
19 // QA objects are 1 & 2 Dimensional histograms. //
20 // author: S.Arcelli //
22 ///////////////////////////////////////////////////////////////////////
25 Modified by fbellini on 17/10/2011
26 - fix for memory leak in constructor
27 - added methods to read histos ranges from config file in DQM
28 - added CTTM maps + relative methods to retrieve CTTM numbering
29 - removed obslete comments
31 Modified by fbellini & rshanoian on 06/07/2011
32 - changes for trigger classes implementation
33 - fRunNumber added as private member
34 - added time vs BCID plot
36 Modified by fbellini on 18/01/2011
37 - reduced histo binning to reduce size
38 - added decoding errors plot
39 - added channel maps and options for DQM shifters
40 - new list of recPoints and ESDs plots
41 - removed hTrmChannels035 and hTrmChannels3671
43 Modified by bvonhall on 03/11/2010
44 - modified declaration of hTrmChannels035 and hTrmChannels3671 in EndOfDetectorCycle()
45 to prevent memory corruption
47 Modified by adecaro on 18/10/2010
48 - fTOFRawStream object set as private member
50 Modified by fbellini on 13/09/2010
51 - Set TLines as private members
52 - Set image flag for expert histos
54 Modified by fbellini on 14/06/2010
56 - use LoadRawDataBuffersV2()
58 Modified by fbellini on 10/05/2010
59 - Fixed EndOfDetectorCycle() memory corruption bug
61 Modified by fbellini on 22/04/2010
62 - Added filter for physics events
64 Modified by fbellini on 16/04/2010
65 - Added EnableDqmShifterOpt()
66 - Modified EndOfDetectorCycle() with options for DQM
69 Modified by fbellini on 30/03/2010
70 - Changed raws time histos range to 610ns
71 - Added FilterLTMData() and FilterSpare() methods
72 - Added check on enabled channels for raw data
73 - Updated RecPoints QA
75 Modified by fbellini on 02/03/2010
76 - Fixed raw data decoding methods (use AliTOFRawStream::LoadRawDataBuffer())
77 - Added filter for noisy channels and read map from OCDB
78 - Added GetCalibData() method
79 - Added CheckVolumeID() and CheckEquipID() methods
83 #include <TClonesArray.h>
87 #include <TPaveText.h>
89 #include "AliCDBManager.h"
90 #include "AliCDBEntry.h"
91 #include "AliESDEvent.h"
92 #include "AliESDtrack.h"
93 #include "AliQAChecker.h"
94 #include "AliRawReader.h"
95 #include "AliTOFRawStream.h"
96 #include "AliTOFcluster.h"
97 #include "AliTOFQADataMakerRec.h"
98 #include "AliTOFrawData.h"
99 #include "AliTOFGeometry.h"
100 #include "AliTOFChannelOnlineStatusArray.h"
101 #include "AliTOFDecoderSummaryData.h"
102 #include "AliTOFDecoderV2.h"
104 ClassImp(AliTOFQADataMakerRec)
106 Int_t AliTOFQADataMakerRec::fgNbinsMultiplicity=200; //number of bins in multiplicity plot
107 Int_t AliTOFQADataMakerRec::fgRangeMinMultiplicity=0;//min range in multiplicity plot
108 Int_t AliTOFQADataMakerRec::fgRangeMaxMultiplicity=200;//max range in multiplicity plot
109 Int_t AliTOFQADataMakerRec::fgNbinsTime=250;//number of bins in time plot
110 const Float_t AliTOFQADataMakerRec::fgkNbinsWidthTime=2.44;//width of bins in time plot
111 Float_t AliTOFQADataMakerRec::fgRangeMinTime=0.0;//range min in time plot
112 Float_t AliTOFQADataMakerRec::fgRangeMaxTime=620.0; //range max in time plot
113 Int_t AliTOFQADataMakerRec::fgCutNmaxFiredMacropad=5;//cut on number of max fired macropad
114 const Int_t AliTOFQADataMakerRec::fgkFiredMacropadLimit=50;//cut on number of max fired macropad
117 //____________________________________________________________________________
118 AliTOFQADataMakerRec::AliTOFQADataMakerRec() :
119 AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kTOF), "TOF Quality Assurance Data Maker"),
121 fEnableNoiseFiltering(kFALSE),
122 fEnableDqmShifterOpt(kFALSE),
124 fLineExpTimeMin(0x0),
125 fLineExpTimeMax(0x0),
128 fTOFRawStream(AliTOFRawStream()),
131 fCalib(AliTOFcalib())
136 // fLineExpTimeMin = new TLine(200., 0., 200., 0.);
137 // fLineExpTimeMax = new TLine(250., 0., 250., 0.);
138 // fLineExpTotMin = new TLine(5., 0., 5., 0.);
139 // fLineExpTotMax = new TLine(20., 0., 20., 0.);
140 for (Int_t sm=0;sm<17;sm++){
141 fLineSMid[sm] = new TLine( sm+1, 0., sm+1, 91.);
144 for (Int_t sm=0;sm<71;sm++){
145 fLineLTMid[sm] = new TLine( sm+1, 0., sm+1, 23.);
148 for (Int_t sm=0;sm<22;sm++){
149 fLineLTMbitId[sm] = new TLine( 0., sm+1, 72. ,sm+1);
154 //____________________________________________________________________________
155 AliTOFQADataMakerRec::AliTOFQADataMakerRec(const AliTOFQADataMakerRec& qadm) :
157 fCalibData(qadm.fCalibData),
158 fEnableNoiseFiltering(qadm.fEnableNoiseFiltering),
159 fEnableDqmShifterOpt(qadm.fEnableDqmShifterOpt),
161 fLineExpTimeMin(qadm.fLineExpTimeMin),
162 fLineExpTimeMax(qadm.fLineExpTimeMax),
163 fLineExpTotMin(qadm.fLineExpTotMin),
164 fLineExpTotMax(qadm.fLineExpTotMax),
165 fTOFRawStream(qadm.fTOFRawStream),
166 fDecoderSummary(qadm.fDecoderSummary),
167 fRunNumber(qadm.fRunNumber),
173 SetName((const char*)qadm.GetName()) ;
174 SetTitle((const char*)qadm.GetTitle());
176 for (Int_t sm=0;sm<17;sm++){
177 fLineSMid[sm]=qadm.fLineSMid[sm];
180 for (Int_t sm=0;sm<71;sm++){
181 fLineLTMid[sm] = qadm.fLineLTMid[sm];
184 for (Int_t sm=0;sm<22;sm++){
185 fLineLTMbitId[sm] = qadm.fLineLTMbitId[sm];
189 //__________________________________________________________________
190 AliTOFQADataMakerRec& AliTOFQADataMakerRec::operator = (const AliTOFQADataMakerRec& qadm )
193 // assignment operator.
195 this->~AliTOFQADataMakerRec();
196 new(this) AliTOFQADataMakerRec(qadm);
200 //----------------------------------------------------------------------------
201 AliTOFQADataMakerRec::~AliTOFQADataMakerRec()
204 fTOFRawStream.Clear();
207 delete fLineExpTimeMin;
209 delete fLineExpTimeMax;
211 delete fLineExpTotMin;
213 delete fLineExpTotMax;
214 for (Int_t sm=0;sm<17;sm++){
216 delete fLineSMid[sm];
218 for (Int_t sm=0;sm<71;sm++){
220 delete fLineLTMid[sm];
222 for (Int_t sm=0;sm<22;sm++){
223 if (fLineLTMbitId[sm])
224 delete fLineLTMbitId[sm];
227 //----------------------------------------------------------------------------
228 AliTOFChannelOnlineStatusArray* AliTOFQADataMakerRec::GetCalibData()
231 // Retrive TOF calib objects from OCDB
233 AliCDBManager *man = AliCDBManager::Instance();
236 if (fRun<=0) fRunNumber=145288; //reference run from LHC11a
237 else fRunNumber=fRun;
239 if (man->GetRun()!=fRunNumber){
240 fRunNumber=man->GetRun();
241 AliWarning(Form("Run number mismatch found: setting it to value from current AliCDBManager instance = %i", fRunNumber));
243 cdbe = man->Get("TOF/Calib/Status",fRunNumber);
247 AliWarning("Load of calibration data from default (alien://) storage failed!");
248 printf("Calibration data will be loaded from local storage - ok if on DQM station!");
249 man->SetDefaultStorage("local:///local/cdb/");
250 cdbe = man->Get("TOF/Calib/Status",fRun);
253 AliWarning("Load of calibration data from local DQM machine storage failed!");
254 AliWarning("Calibration data will be loaded from local ($ALICE_ROOT) storage ");
255 man->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
256 cdbe = man->Get("TOF/Calib/Status",fRunNumber);
259 // Retrieval of data in directory TOF/Calib/Data:
260 AliTOFChannelOnlineStatusArray * array = 0;
262 printf("======= OCDB object for TOF retrieved from run %i in %s\n",fRunNumber,cdbe->GetName());
263 array = (AliTOFChannelOnlineStatusArray *)cdbe->GetObject();
265 if (!array) AliFatal("No calibration data from calibration database !");
267 fCalib.Init(fRunNumber);
271 //____________________________________________________________________________
272 void AliTOFQADataMakerRec::InitRaws()
275 // create Raws histograms in Raws subdir
277 ReadHistogramRangeFromFile(gSystem->Getenv("TOFDQMHISTO_CONFIGFILE"));
279 const Bool_t expert = kTRUE ;
280 const Bool_t saveCorr = kTRUE ;
281 const Bool_t image = kTRUE ;
283 TH1I * h0 = new TH1I("hTOFRaws","TOF raw hit multiplicity; TOF raw hits number; Events ",fgNbinsMultiplicity, fgRangeMinMultiplicity, fgRangeMaxMultiplicity);
284 TH1I * h1 = new TH1I("hTOFRawsIA","TOF raw hit multiplicity - I/A side; TOF raw hits number;Events ",fgNbinsMultiplicity, fgRangeMinMultiplicity, fgRangeMaxMultiplicity);
285 TH1I * h2 = new TH1I("hTOFRawsOA","TOF raw hit multiplicity - O/A side; TOF raw hits number;Events ",fgNbinsMultiplicity, fgRangeMinMultiplicity, fgRangeMaxMultiplicity);
286 TH1I * h3 = new TH1I("hTOFRawsIC","TOF raw hit multiplicity - I/C side; TOF raw hits number;Events ",fgNbinsMultiplicity, fgRangeMinMultiplicity, fgRangeMaxMultiplicity);
287 TH1I * h4 = new TH1I("hTOFRawsOC","TOF raw hit multiplicity - O/C side; TOF raw hits number;Events ",fgNbinsMultiplicity, fgRangeMinMultiplicity, fgRangeMaxMultiplicity);
289 TH1F * h5 = new TH1F("hTOFRawsTime", "TOF Raws - Hit time (ns);Measured Hit time [ns];Hits", fgNbinsTime,fgRangeMinTime,fgRangeMaxTime);
290 TH1F * h6 = new TH1F("hTOFRawsTimeIA", "TOF Raws - Hit time (ns) - I/A side;Measured Hit time [ns];Hits", fgNbinsTime,fgRangeMinTime,fgRangeMaxTime);
291 TH1F * h7 = new TH1F("hTOFRawsTimeOA", "TOF Raws - Hit time (ns) - O/A side;Measured Hit time [ns];Hits", fgNbinsTime,fgRangeMinTime,fgRangeMaxTime);
292 TH1F * h8 = new TH1F("hTOFRawsTimeIC", "TOF Raws - Hit time (ns) - I/C side;Measured Hit time [ns];Hits", fgNbinsTime,fgRangeMinTime,fgRangeMaxTime);
293 TH1F * h9 = new TH1F("hTOFRawsTimeOC", "TOF Raws - Hit time (ns) - O/C side;Measured Hit time [ns];Hits", fgNbinsTime,fgRangeMinTime,fgRangeMaxTime);
295 TH1F * h10 = new TH1F("hTOFRawsToT", "TOF Raws - Hit ToT (ns);Measured Hit ToT (ns);Hits", 100, 0., 48.8);
297 TH1F * h11 = new TH1F("hTOFRawsToTIA", "TOF Raws - Hit ToT (ns) - I/A side;Measured Hit ToT (ns);Hits", 100, 0., 48.8);
298 TH1F * h12 = new TH1F("hTOFRawsToTOA", "TOF Raws - Hit ToT (ns) - O/A side;Measured Hit ToT (ns);Hits", 100, 0., 48.8);
299 TH1F * h13 = new TH1F("hTOFRawsToTIC", "TOF Raws - Hit ToT (ns) - I/C side;Measured Hit ToT (ns);Hits", 100, 0., 48.8);
300 TH1F * h14 = new TH1F("hTOFRawsToTOC", "TOF Raws - Hit ToT (ns) - O/C side;Measured Hit ToT (ns);Hits", 100, 0., 48.8);
302 TH1F * h15 = new TH1F("hTOFRawsLTMHits", "LTMs OR signals; Crate; Counts", 72, 0., 72.);
303 TH2F * h16 = new TH2F("hTOFrefMap", "TOF enabled channel reference map;sector;strip", 72, 0., 18., 91, 0., 91.);
304 TH2F * h17 = new TH2F("hTOFRawHitMap","TOF raw hit map;sector;strip", 72, 0., 18., 91, 0., 91.);
305 TH2I * h18 = new TH2I("hTOFDecodingErrors","Decoding error monitoring; DDL; Error ", 72, 0, 72, 13,1,14);
307 h18->GetYaxis()->SetBinLabel(1,"DRM ");
308 h18->GetYaxis()->SetBinLabel(2,"LTM ");
309 h18->GetYaxis()->SetBinLabel(3,"TRM 3 ");
310 h18->GetYaxis()->SetBinLabel(4,"TRM 4");
311 h18->GetYaxis()->SetBinLabel(5,"TRM 5");
312 h18->GetYaxis()->SetBinLabel(6,"TRM 6");
313 h18->GetYaxis()->SetBinLabel(7,"TRM 7");
314 h18->GetYaxis()->SetBinLabel(8,"TRM 8");
315 h18->GetYaxis()->SetBinLabel(9,"TRM 9");
316 h18->GetYaxis()->SetBinLabel(10,"TRM 10");
317 h18->GetYaxis()->SetBinLabel(11,"TRM 11");
318 h18->GetYaxis()->SetBinLabel(12,"TRM 12");
319 h18->GetYaxis()->SetBinLabel(13,"recovered");
321 TH1F * h19 = new TH1F("hTOFOrphansTime", "TOF Raws - Orphans time (ns);Measured Hit time [ns];Hits",fgNbinsTime,fgRangeMinTime,fgRangeMaxTime);
322 TH2F * h20 = new TH2F("hTOFRawTimeVsTRM035", "TOF raws - Hit time vs TRM - crates 0 to 35; TRM index = DDL*10+TRM(0-9);TOF raw time [ns]", 361, 0., 361.,fgNbinsTime,fgRangeMinTime,fgRangeMaxTime);
323 TH2F * h21 = new TH2F("hTOFRawTimeVsTRM3671", "TOF raws - Hit time vs TRM - crates 36 to 72; TRM index = DDL**10+TRM(0-9);TOF raw time [ns]", 361, 360., 721.,fgNbinsTime, fgRangeMinTime,fgRangeMaxTime);
324 TH2F * h22 = new TH2F("hTOFTimeVsStrip","TOF raw hit time vs. MRPC (along z axis); MRPC index along z axis; Raws TOF time (ns) ", 91,0.,91,fgNbinsTime,fgRangeMinTime,fgRangeMaxTime);
325 TH2F * h23 = new TH2F("hTOFtimeVsBCID","TOF time vs BCID; BCID; time (ns) ", 3564, 0., 3564.,fgNbinsTime,fgRangeMinTime,fgRangeMaxTime);
326 TH2F * h24 = new TH2F("hTOFchannelEfficiencyMap","TOF channels (HWok && efficient && !noisy && !problematic);sector;strip", 72, 0., 18., 91, 0., 91.);
327 TH2F * h25 = new TH2F("hTOFhitsCTTM","Map of hit pads according to CTTM numbering;LTM index;bit index", 72, 0., 72., 23, 0., 23.);
328 TH2F * h26 = new TH2F("hTOFmacropadCTTM","Map of hit macropads according to CTTM numbering;LTM index; bit index", 72, 0., 72., 23, 0., 23.);
329 TH2F * h27 = new TH2F("hTOFmacropadDeltaPhiTime","#Deltat vs #Delta#Phi of hit macropads;#Delta#Phi (degrees);#DeltaBX", 18, 0., 180., 20, 0., 20.0);
331 h25->GetYaxis()->SetTickLength(-0.02);
332 h26->GetYaxis()->SetTickLength(-0.02);
333 h25->GetYaxis()->SetNdivisions(210);
334 h26->GetYaxis()->SetNdivisions(210);
335 h25->GetXaxis()->SetTickLength(-0.02);
336 h26->GetXaxis()->SetTickLength(-0.02);
337 h25->GetXaxis()->SetLabelOffset(0.015);
338 h26->GetXaxis()->SetLabelOffset(0.015);
339 h25->GetXaxis()->SetNdivisions(515);
340 h26->GetXaxis()->SetNdivisions(515);
370 //add lines for DQM shifter
371 fLineExpTimeMin = new TLine(200., 0., 200., 0.);
372 fLineExpTimeMax = new TLine(250., 0., 250., 0.);
373 fLineExpTotMin = new TLine(5., 0., 5., 0.);
374 fLineExpTotMax = new TLine(20., 0., 20., 0.);
376 fLineExpTimeMin->SetLineColor(kGreen);
377 fLineExpTimeMin->SetLineWidth(2);
379 fLineExpTimeMax->SetLineColor(kGreen);
380 fLineExpTimeMax->SetLineWidth(2);
382 fLineExpTotMin->SetLineColor(kGreen);
383 fLineExpTotMin->SetLineWidth(2);
385 fLineExpTotMax->SetLineColor(kGreen);
386 fLineExpTotMax->SetLineWidth(2);
388 for (Int_t sm=0;sm<17;sm++){
389 fLineSMid[sm]->SetLineColor(kMagenta);
390 fLineSMid[sm]->SetLineWidth(2);
393 h5->GetListOfFunctions()->Add(fLineExpTimeMin);
394 h5->GetListOfFunctions()->Add(fLineExpTimeMax);
395 h10->GetListOfFunctions()->Add(fLineExpTotMin);
396 h10->GetListOfFunctions()->Add(fLineExpTotMax);
398 for (Int_t sm=0;sm<17;sm++){
399 h16->GetListOfFunctions()->Add(fLineSMid[sm]);
400 h17->GetListOfFunctions()->Add(fLineSMid[sm]);
403 for (Int_t sm=0;sm<71;sm++){
404 fLineLTMid[sm]->SetLineColor(kBlack);
405 fLineLTMid[sm]->SetLineWidth(1);
406 h26->GetListOfFunctions()->Add(fLineLTMid[sm]);
407 h25->GetListOfFunctions()->Add(fLineLTMid[sm]);
409 for (Int_t sm=0;sm<22;sm++){
410 fLineLTMbitId[sm]->SetLineColor(kBlack);
411 fLineLTMbitId[sm]->SetLineWidth(1);
412 h26->GetListOfFunctions()->Add(fLineLTMbitId[sm]);
413 h25->GetListOfFunctions()->Add(fLineLTMbitId[sm]);
416 TPaveText *phosHoleBox=new TPaveText(13,38,16,53,"b");
417 phosHoleBox->SetFillStyle(0);
418 phosHoleBox->SetFillColor(kWhite);
419 phosHoleBox->SetLineColor(kMagenta);
420 phosHoleBox->SetLineWidth(2);
421 phosHoleBox->AddText("PHOS");
422 h16->GetListOfFunctions()->Add(phosHoleBox);
423 h17->GetListOfFunctions()->Add(phosHoleBox);
425 // h0->SetDrawOption("logy");
426 // h5->SetDrawOption("logy");
427 // h10->SetDrawOption("logy");
429 Add2RawsList(h0, 0, !expert, image, !saveCorr) ;
430 Add2RawsList(h1, 1, expert, !image, !saveCorr) ;
431 Add2RawsList(h2, 2, expert, !image, !saveCorr) ;
432 Add2RawsList(h3, 3, expert, !image, !saveCorr) ;
433 Add2RawsList(h4, 4, expert, !image, !saveCorr) ;
434 Add2RawsList(h5, 5, !expert, image, !saveCorr) ;
435 Add2RawsList(h6, 6, expert, !image, !saveCorr) ;
436 Add2RawsList(h7, 7, expert, !image, !saveCorr) ;
437 Add2RawsList(h8, 8, expert, !image, !saveCorr) ;
438 Add2RawsList(h9, 9, expert, !image, !saveCorr) ;
439 Add2RawsList(h10, 10, !expert, image, !saveCorr) ;
440 Add2RawsList(h11, 11, expert, !image, !saveCorr) ;
441 Add2RawsList(h12, 12, expert, !image, !saveCorr) ;
442 Add2RawsList(h13, 13, expert, !image, !saveCorr) ;
443 Add2RawsList(h14, 14, expert, !image, !saveCorr) ;
444 Add2RawsList(h15, 15, expert, !image, !saveCorr) ;
445 Add2RawsList(h16, 16, !expert, image, !saveCorr) ;
446 Add2RawsList(h17, 17, !expert, image, !saveCorr) ;
447 Add2RawsList(h18, 18, expert, !image, !saveCorr) ;
448 Add2RawsList(h19, 19, expert, !image, !saveCorr) ;
449 Add2RawsList(h20, 20, expert, !image, !saveCorr) ;
450 Add2RawsList(h21, 21, expert, !image, !saveCorr) ;
451 Add2RawsList(h22, 22, !expert, image, !saveCorr) ;
452 Add2RawsList(h23, 23, !expert, !image, !saveCorr) ;
453 Add2RawsList(h24, 24, !expert, !image, !saveCorr) ;
454 Add2RawsList(h25, 25, !expert, !image, !saveCorr) ;
455 Add2RawsList(h26, 26, !expert, image, !saveCorr) ;
456 Add2RawsList(h27, 27, !expert, image, !saveCorr) ;
458 ClonePerTrigClass(AliQAv1::kRAWS); // this should be the last line
461 //____________________________________________________________________________
462 void AliTOFQADataMakerRec::InitRecPoints()
465 // create RecPoints histograms in RecPoints subdir
468 const Bool_t expert = kTRUE ;
469 const Bool_t image = kTRUE ;
471 TH1I * h0 = new TH1I("hTOFRecPoints", "TOF RecPoints multiplicity ; TOF RecPoints number;Events",200, 0, 200) ;
473 TH1F * h1 = new TH1F("hTOFRecPointsTimeIA", "RecPoints Time Spectrum in TOF (ns)- I/A side; Calibrated TOF time [ns];Events", 250, 0., 610.) ;
474 TH1F * h2 = new TH1F("hTOFRecPointsTimeOA", "RecPoints Time Spectrum in TOF (ns)- O/A side; Calibrated TOF time [ns];Events", 250, 0., 610.) ;
475 TH1F * h3 = new TH1F("hTOFRecPointsTimeIC", "RecPoints Time Spectrum in TOF (ns)- I/C side; Calibrated TOF time [ns];Events", 250, 0., 610.) ;
476 TH1F * h4 = new TH1F("hTOFRecPointsTimeOC", "RecPoints Time Spectrum in TOF (ns)- O/C side; Calibrated TOF time [ns];Events", 250, 0., 610.) ;
478 TH1F * h5 = new TH1F("hTOFRecPointsRawTimeIA", "RecPoints raw Time Spectrum in TOF (ns)-I/A side; Measured TOF time [ns];Hits", 250, 0., 610.) ;
479 TH1F * h6 = new TH1F("hTOFRecPointsRawTimeOA", "RecPoints raw Time Spectrum in TOF (ns)-O/A side; Measured TOF time [ns];Hits", 250, 0., 610.) ;
480 TH1F * h7 = new TH1F("hTOFRecPointsRawTimeIC", "RecPoints raw Time Spectrum in TOF (ns)-I/C side; Measured TOF time [ns];Hits", 250, 0., 610.) ;
481 TH1F * h8 = new TH1F("hTOFRecPointsRawTimeOC", "RecPoints raw Time Spectrum in TOF (ns)-O/C side; Measured TOF time [ns];Hits", 250, 0., 610.) ;
483 TH1F * h9 = new TH1F("hTOFRecPointsToTIA", "RecPoints ToT Spectrum in TOF (ns)-I/A side; Measured TOT [ns];Hits", 100, 0., 48.8 ) ;
484 TH1F * h10 = new TH1F("hTOFRecPointsToTOA", "RecPoints ToT Spectrum in TOF (ns)-O/A side; Measured TOT [ns];Hits", 100, 0., 48.8 ) ;
485 TH1F * h11 = new TH1F("hTOFRecPointsToTIC", "RecPoints ToT Spectrum in TOF (ns)-I/C side; Measured TOT [ns];Hits", 100, 0., 48.8 ) ;
486 TH1F * h12 = new TH1F("hTOFRecPointsToTOC", "RecPoints ToT Spectrum in TOF (ns)-O/C side; Measured TOT [ns];Hits", 100, 0., 48.8 ) ;
488 TH2F * h13 = new TH2F("hTOFRecPointsClusMap","RecPoints map; sector ;strip", 72, 0., 18., 91, 0., 91.) ;
489 TH2F * h14 = new TH2F("hTOFRecPointsTimeVsStrip","RecPoints TOF time vs. strip (theta); Strip index; RecPoints TOF time (ns) ",91, 0., 91., 250, 0., 610.) ;
490 TH2F * h15 = new TH2F("hTOFRecPointsTimeVsTRM035","TOF RecPoints time vs TRM - crates 0 to 35; TRM index = DDL*10+TRM(0-9);TOF time [ns]", 361, 0., 361., 250, 0., 610.0) ;
491 TH2F * h16 = new TH2F("hTOFRecPointsTimeVsTRM3671","TOF RecPoints time vs TRM - crates 36 to 72; TRM index = DDL**10+TRM(0-9);TOF time [ns]", 361, 360., 721., 250, 0., 610.0) ;
511 Add2RecPointsList(h0, 0, !expert, image) ;
512 Add2RecPointsList(h1, 1, !expert, image) ;
513 Add2RecPointsList(h2, 2, !expert, image) ;
514 Add2RecPointsList(h3, 3, !expert, image) ;
515 Add2RecPointsList(h4, 4, !expert, image) ;
516 Add2RecPointsList(h5, 5, expert, !image) ;
517 Add2RecPointsList(h6, 6, expert, !image) ;
518 Add2RecPointsList(h7, 7, expert, !image) ;
519 Add2RecPointsList(h8, 8, expert, !image) ;
520 Add2RecPointsList(h9, 9, !expert, !image) ;
521 Add2RecPointsList(h10, 10, !expert, !image) ;
522 Add2RecPointsList(h11, 11, !expert, !image) ;
523 Add2RecPointsList(h12, 12, !expert, !image) ;
524 Add2RecPointsList(h13, 13, expert, !image) ;
525 Add2RecPointsList(h14, 14, expert, image) ;
526 Add2RecPointsList(h15, 15, expert, !image) ;
527 Add2RecPointsList(h16, 16, expert, !image) ;
529 ClonePerTrigClass(AliQAv1::kRECPOINTS); // this should be the last line
532 //____________________________________________________________________________
533 void AliTOFQADataMakerRec::InitESDs()
536 //create ESDs histograms in ESDs subdir
539 const Bool_t expert = kTRUE ;
540 const Bool_t image = kTRUE ;
542 TH1I * h0 = new TH1I("hTOFESDs", "Number of matched TOF tracks per event;Number of TOF matched ESD tracks;Counts", 200, 0, 200) ;
543 TH1F * h1 = new TH1F("hTOFESDsTime", "Matched ESDs tracks: TOF Time spectrum; Calibrated TOF time [ns];Counts", 250, 0., 610. ) ;
544 TH1F * h2 = new TH1F("hTOFESDsRawTime", "Matched ESDs tracks: TOF raw Time spectrum;Measured TOF time [ns];Counts", 250, 0., 610.) ;
545 TH1F * h3 = new TH1F("hTOFESDsToT", "Matched ESDs tracks: TOF ToT spectrum; ESDs ToT [ns];Counts",100, 0., 48.8) ;
546 TH1F * h4 = new TH1F("hTOFESDskTOFOUT", "p_{T} distribution of tracks with kTOFout; p_{T} (GeV/c);Counts", 50, 0.20, 5.00) ;
547 TH1F * h5 = new TH1F("hTOFESDskTIME", "p_{T} distribution of tracks with kTOFout && kTIME; p_{T} (GeV/c);Counts", 50, 0.20, 5.00) ;
548 TH1F * h6 = new TH1F("hTOFESDsMatched", "p_{T} distribution of tracks with kTOFout && TOFtime>0; p_{T} (GeV/c);Counts", 50, 0.20, 5.00) ;
549 TH1F * h7 = new TH1F("hTOFESDsMatchingProb", "TPC-TOF track-matching probability;TOF matching probability (%) ;Counts",101, -1.0, 100) ;
550 TH1F * h8 = new TH1F("hTOFESDsDiffTime", "ESDs t_{TOF}-t_{exp,pi} spectrum in TOF (ps); t_{TOF}-t_{exp,pi} [ps];Counts", 200, -2440., 2440.) ;
551 TH1F * h9 = new TH1F("hTOFHitsLength", "Matched ESDs tracks: Length Spectrum; Track length [cm];Counts", 800, 0., 800) ;
564 Add2ESDsList(h0, 0, !expert, image) ;
565 Add2ESDsList(h1, 1, !expert, image) ;
566 Add2ESDsList(h2, 2, expert, !image) ;
567 Add2ESDsList(h3, 3, !expert, !image) ;
568 Add2ESDsList(h4, 4, expert, image) ;
569 Add2ESDsList(h5, 5, expert, image) ;
570 Add2ESDsList(h6, 6, expert, image) ;
571 Add2ESDsList(h7, 7, expert, image) ;
572 Add2ESDsList(h8, 8, expert, !image) ;
573 Add2ESDsList(h9, 9, !expert, !image) ;
575 ClonePerTrigClass(AliQAv1::kESDS); // this should be the last line
579 //____________________________________________________________________________
580 void AliTOFQADataMakerRec::MakeRaws(AliRawReader* rawReader)
583 // makes data from Raws
585 if (rawReader->GetType()==7) {
587 Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;//in ns
588 Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
589 Int_t ntof[5]; /* 0=tot, 1=IA, 2=OA, 3=IC, 4=OC*/
590 for (Int_t j=0;j<5;j++){ ntof[j]=0;}
591 Int_t equipmentID[5]; //(ddl, trm, chain,tdc,channel)
592 Int_t volumeID[5]; //(sector,plate,strip,padX,padZ)
593 Int_t volumeID2[5]; //(sector,plate,strip,padZ,padX) to use AliTOFGeometry::GetIndex()
595 Int_t indexCTTM[2]={-1,-1};
596 Int_t indexGeo2CTTM[2]={-1,-1};
597 Float_t macropadPhiTimeUPC[fgkFiredMacropadLimit][2];
598 for (Int_t ii=0;ii<2;ii++){
599 for (Int_t jj=0;jj<fgkFiredMacropadLimit;jj++){
600 macropadPhiTimeUPC[jj][ii]=-999.0;
604 TClonesArray * clonesRawData;
605 fTOFRawStream.SetRawReader(rawReader);
606 Int_t BCID=rawReader->GetBCID();
608 Int_t nFiredMacropad=0,
610 nFiredMacropad=GetNumberOfFiredMacropad(rawReader);
612 //uncomment if needed to apply DeltaBC correction
613 //fTOFRawStream.ApplyBCCorrections(kTRUE);
615 if (fDecoderSummary){
616 fDecoderSummary->Reset();
618 for (Int_t iDDL = 0; iDDL < AliTOFGeometry::NDDL()*AliTOFGeometry::NSectors(); iDDL++){
620 fTOFRawStream.LoadRawDataBuffersV2(iDDL);
622 //* get decoding error counters
623 fDecoderSummary = ( (AliTOFDecoderV2*) fTOFRawStream.GetDecoderV2() )->GetDecoderSummaryData();
624 if ( (fDecoderSummary) && (fDecoderSummary ->GetErrorDetected()) ) {
625 Int_t errorSlotID=(Int_t) fDecoderSummary->GetErrorSlotID();
626 FillRawsData(18,iDDL,errorSlotID);
627 if (fDecoderSummary -> GetRecoverError() )
628 FillRawsData(18,iDDL,13);
631 clonesRawData = (TClonesArray*)fTOFRawStream.GetRawData();
632 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
633 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
635 if (tofRawDatum->GetTOF()){
637 equipmentID[1]=tofRawDatum->GetTRM();
638 equipmentID[2]=tofRawDatum->GetTRMchain();
639 equipmentID[3]=tofRawDatum->GetTDC();
640 equipmentID[4]=tofRawDatum->GetTDCchannel();
642 if (CheckEquipID(equipmentID)){
643 fTOFRawStream.EquipmentId2VolumeId(iDDL,
644 tofRawDatum->GetTRM(),
645 tofRawDatum->GetTRMchain(),
646 tofRawDatum->GetTDC(),
647 tofRawDatum->GetTDCchannel(),
650 if (FilterLTMData(equipmentID)) { //counts LTM hits
651 if (equipmentID[2]==1) //crate left, A-side or C-side
652 FillRawsData(15,equipmentID[0]);
654 FillRawsData(15,equipmentID[0]-1);
656 //fired macropad map (from LTM hits) - only for low multi evts (UPC)
657 if ((nFiredMacropad<=fgCutNmaxFiredMacropad)){
659 AliInfo(Form("Event found with %i fired macropads in BCID = %i!",nFiredMacropad,BCID));
661 GetCTTMIndex(equipmentID, indexCTTM);
662 FillRawsData(26,indexCTTM[0],indexCTTM[1]);
663 Float_t halfSMphi=-999.0;
666 halfSMphi=indexCTTM[0]*10.+5.;
667 else halfSMphi=(indexCTTM[0]-36)*10.+5.;
668 macropadPhiTimeUPC[iFiredMacropad][0]=halfSMphi;
669 indexBC= TMath::Nint(tofRawDatum->GetTOF()*tdc2ns)/25;
670 macropadPhiTimeUPC[iFiredMacropad][1]=indexBC;
675 if (CheckVolumeID(volumeID)){
676 volumeID2[0]=volumeID[0];
677 volumeID2[1]=volumeID[1];
678 volumeID2[2]=volumeID[2];
679 volumeID2[3]=volumeID[4];
680 volumeID2[4]=volumeID[3];
681 chIndex=AliTOFGeometry::GetIndex(volumeID2);
683 //fill hit map according to CTTM numbering
684 GetGeo2CTTMIndex(volumeID2, indexGeo2CTTM);
685 if ((nFiredMacropad<=fgCutNmaxFiredMacropad)){
686 FillRawsData(25,indexGeo2CTTM[0],indexGeo2CTTM[1]);
689 if (tofRawDatum->GetTOT()){
690 if (!(fCalibData->GetNoiseStatus(chIndex)==AliTOFChannelOnlineStatusArray::kTOFNoiseBad)
691 && (fCalibData->GetHWStatus(chIndex) == AliTOFChannelOnlineStatusArray::kTOFHWOk)) {//noise and enabled filter
692 ntof[0]++; //counter for tof hits
694 //fill global spectra for DQM plots
695 FillRawsData(5, tofRawDatum->GetTOF()*tdc2ns) ;//in ns
696 FillRawsData(10, tofRawDatum->GetTOT()*tot2ns) ;//in ns
697 FillRawsData(23, BCID, tofRawDatum->GetTOF()*tdc2ns) ;//in ns
699 //fill side-related spectra for experts plots
700 Int_t ddlACside=iDDL/36; // 0 or 1
701 Int_t ddlPerSm=iDDL%4;
703 if (volumeID2[0]>4 && volumeID2[0]<14){ //O side
704 if (ddlPerSm<2){ //A side
706 FillRawsData(6, tofRawDatum->GetTOF()*tdc2ns) ;
707 FillRawsData(11, tofRawDatum->GetTOT()*tot2ns) ;
710 FillRawsData(8, tofRawDatum->GetTOF()*tdc2ns) ;
711 FillRawsData(13, tofRawDatum->GetTOT()*tot2ns) ;
714 if (volumeID2[0]<5 || volumeID2[0]>13){ //I side
715 if (ddlPerSm<2){ //A side
717 FillRawsData(7, tofRawDatum->GetTOF()*tdc2ns) ;
718 FillRawsData(12, tofRawDatum->GetTOT()*tot2ns) ;
721 FillRawsData(9, tofRawDatum->GetTOF()*tdc2ns) ;
722 FillRawsData(14, tofRawDatum->GetTOT()*tot2ns) ;
728 Int_t trm= iDDL*10+(equipmentID[1]-3);
729 FillRawsData(20+ddlACside,trm,tofRawDatum->GetTOF()*tdc2ns);
730 FillRawsData(22,GetStripIndex(volumeID),tofRawDatum->GetTOF()*tdc2ns) ;
731 Short_t fea = volumeID2[4]/12;
732 Float_t hitmapx = volumeID2[0] + ((Double_t)(3 - fea) + 0.5) *0.25;
733 FillRawsData(17,hitmapx,GetStripIndex(volumeID2));
734 if (fCalib.IsChannelEnabled(chIndex,kTRUE,kTRUE))//checks also if efficient and if problematic
735 FillRawsData(24,hitmapx,GetStripIndex(volumeID2));
739 if (!(fCalibData->GetNoiseStatus(chIndex) == AliTOFChannelOnlineStatusArray::kTOFNoiseBad)
740 && (fCalibData->GetHWStatus(chIndex) == AliTOFChannelOnlineStatusArray::kTOFHWOk))
741 FillRawsData(19, tofRawDatum->GetTOF()*tdc2ns) ;//in ns
743 }//end volumeID check
747 clonesRawData->Clear();
750 for (Int_t j=0;j<5;j++) FillRawsData(j,ntof[j]);
751 fTOFRawStream.Clear();
753 if ((nFiredMacropad<=fgCutNmaxFiredMacropad)){
754 Float_t deltaPhiMacropad=-999.;
755 Float_t deltaTimeMacropad=-999.;
756 for (Int_t j=0;j<fgCutNmaxFiredMacropad+1; j++){
757 for (Int_t k=j+1;k<fgCutNmaxFiredMacropad+1; k++){
758 if ((macropadPhiTimeUPC[j][0]>0.0)&&(macropadPhiTimeUPC[k][0]>0.0)){
759 deltaPhiMacropad=TMath::Abs(macropadPhiTimeUPC[j][0]-macropadPhiTimeUPC[k][0]);
760 deltaTimeMacropad=TMath::Abs(macropadPhiTimeUPC[j][1]-macropadPhiTimeUPC[k][1]);
761 if (deltaPhiMacropad<=180.)
762 FillRawsData(27, deltaPhiMacropad,deltaTimeMacropad);
764 FillRawsData(27, TMath::Abs(360.0-deltaPhiMacropad),deltaTimeMacropad);
768 }//end cut on number of fired macropad
770 AliDebug(1,Form("Event of type %d found. Skipping non-physics event for QA.\n", rawReader->GetType()));
773 //fill reference map for DQM shifter only once in a detector cycle
775 Int_t geoId[5]={-1,-1,-1,-1,-1};// pgeoId=(sm, mod, strip, padZ, padX)
776 Int_t detId[5]={-1,-1,-1,-1,-1};//detID=(ddl,trm,tdc, chain,channel)
778 for (Int_t ch = 0; ch < fCalibData->GetSize(); ch++) {
779 AliTOFGeometry::GetVolumeIndices(ch,geoId);
780 AliTOFRawStream::Geant2EquipmentId(geoId,detId);
781 if ((detId[1]<0)||(detId[0]<0)) continue;
782 trmIndex=(detId[1]-3)+detId[0]*10;
784 if ( (!(fCalibData->GetNoiseStatus(ch)==AliTOFChannelOnlineStatusArray::kTOFNoiseBad))
785 && (fCalibData->GetHWStatus(ch) == AliTOFChannelOnlineStatusArray::kTOFHWOk) ){
786 //fill reference map with info from OCDB
787 Short_t fea = geoId[4]/12;
788 Float_t hitmapx = geoId[0] + ((Double_t)(3 - fea) + 0.5)*0.25;
789 FillRawsData(16,hitmapx, GetStripIndex(geoId));
795 //enable options for DQM shifter
796 EnableDqmShifterOpt(kTRUE);
798 IncEvCountCycleRaws();
799 IncEvCountTotalRaws();
803 //____________________________________________________________________________
804 void AliTOFQADataMakerRec::MakeRecPoints(TTree * clustersTree)
807 // Make data from Clusters
810 Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
811 Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
813 Int_t volumeID2[5];//(sm, plate,strip, padZ,padX)
814 //Int_t volumeID[5];//(sm, plate,strip, padX,padZ)
816 TBranch *branch=clustersTree->GetBranch("TOF");
818 AliError("can't get the branch with the TOF clusters !");
822 static TClonesArray dummy("AliTOFcluster",10000);
824 TClonesArray *clusters=&dummy;
825 branch->SetAddress(&clusters);
828 clustersTree->GetEvent(0);
830 FillRecPointsData(0,(Int_t)clusters->GetEntriesFast()) ;
832 TIter next(clusters) ;
834 while ( (c = dynamic_cast<AliTOFcluster *>(next())) ) {
836 // volumeID2[0] = c->GetDetInd(0);
837 // volumeID2[1] = c->GetDetInd(1);
838 // volumeID2[2] = c->GetDetInd(2);
839 // volumeID2[3] = c->GetDetInd(4); //padX
840 // volumeID2[4] = c->GetDetInd(3); //padZ
842 for (Int_t i=0;i<5;i++){
843 volumeID2[i]=c->GetDetInd(i); //X and Z indeces inverted in RecPoints
845 //Int_t chIndex=AliTOFGeometry::GetIndex(volumeID2);
846 Int_t iDDL=AliTOFRawStream::Geant2DDL(volumeID2);
847 Int_t iTRM=AliTOFRawStream::Geant2TRM(volumeID2);
848 Short_t fea = volumeID2[4]/12;
849 Float_t hitmapx = volumeID2[0] + ((Double_t)(3 - fea) + 0.5) *0.25;
850 Int_t ddlACside=iDDL/36; // 0 or 1
851 Int_t ddlPerSm=iDDL%4;
853 if ((c->GetTDCRAW()) && (c->GetTDC()) && (c->GetToT())){
854 if (volumeID2[0]>4 && volumeID2[0]<14){ //I side
855 if (ddlPerSm<2){ //A side
856 FillRecPointsData(1, c->GetTDC()*tdc2ns) ;//in ns
857 FillRecPointsData(5, c->GetTDCRAW()*tdc2ns) ;//in ns
858 FillRecPointsData(9, c->GetToT()*tot2ns) ;//in ns
860 FillRecPointsData(3, c->GetTDC()*tdc2ns) ;//in ns
861 FillRecPointsData(7, c->GetTDCRAW()*tdc2ns) ;//in ns
862 FillRecPointsData(11, c->GetToT()*tot2ns) ;//in ns
865 if (volumeID2[0]<5 || volumeID2[0]>13){ //O side
866 if (ddlPerSm<2){ //A side
867 FillRecPointsData(2, c->GetTDC()*tdc2ns) ;//in ns
868 FillRecPointsData(6, c->GetTDCRAW()*tdc2ns) ;//in ns
869 FillRecPointsData(10, c->GetToT()*tot2ns) ;//in ns
871 FillRecPointsData(4, c->GetTDC()*tdc2ns) ;//in ns
872 FillRecPointsData(8, c->GetTDCRAW()*tdc2ns) ;//in ns
873 FillRecPointsData(12, c->GetToT()*tot2ns) ;//in ns
877 FillRecPointsData(13,hitmapx,GetStripIndex(volumeID2));
878 FillRecPointsData(14,GetStripIndex(volumeID2), c->GetTDC()*tdc2ns) ;
879 Int_t trm= iDDL*10+(iTRM-3);
880 if (ddlACside==0) { //A side
881 FillRecPointsData(15,trm,c->GetTDC()*tdc2ns);
883 FillRecPointsData(16,trm,c->GetTDC()*tdc2ns);
887 EnableDqmShifterOpt(kFALSE);
889 IncEvCountCycleRecPoints();
890 IncEvCountTotalRecPoints();
894 //____________________________________________________________________________
895 void AliTOFQADataMakerRec::MakeESDs(AliESDEvent * esd)
898 // make QA data from ESDs
900 const Double_t speedOfLight = TMath::C()*1E2*1E-12; // cm/ps
901 const Double_t pionMass = 0.13957018; //GeV/c^2
903 Int_t ntrk = esd->GetNumberOfTracks() ;
908 AliESDtrack *track=esd->GetTrack(ntrk);
909 Double_t tofTime=track->GetTOFsignal();//in ps
910 Double_t tofTimeRaw=track->GetTOFsignalRaw();//in ps
911 Double_t tofToT=track->GetTOFsignalToT(); //in ps
913 UInt_t status=track->GetStatus();
914 if (track->IsOn(AliESDtrack::kTPCrefit)) {
916 Double_t y=track->Eta();
917 if (TMath::Abs(y)<0.9) { //select TOF acceptance
918 if ((status&AliESDtrack::kTOFout)!=0) { //define matching
920 FillESDsData(1,tofTime*1E-3);
921 FillESDsData(2,tofTimeRaw*1E-3);
922 FillESDsData(3,tofToT*1E-3);
923 FillESDsData(4,track->Pt());
925 Double_t length =track->GetIntegratedLength();
926 Double_t mom2=(track->Pt()*track->Pt())+(track->Pz()*track->Pz());
927 Double_t piTexp = TMath::Sqrt(1+(pionMass*pionMass/mom2))*length/speedOfLight; //in ps
928 FillESDsData(8,tofTime-piTexp);
929 FillESDsData(9,length);
931 if ((status&AliESDtrack::kTIME)!=0)
932 FillESDsData(5,track->Pt());
935 FillESDsData(6,track->Pt());
936 } //end check on matched tracks
938 }//end check on TPCrefit
941 FillESDsData(0,ntofout) ;
943 Float_t ratio = (Float_t)ntofout/(Float_t)ntpc*100.; //matching probability
944 FillESDsData(7,ratio) ;
948 Float_t ratio = (Float_t)ntofout/(Float_t)ntpc*100; //matched over propagated to TOF outer radius
949 FillESDsData(8,ratio) ;
951 EnableDqmShifterOpt(kFALSE);
953 IncEvCountCycleESDs();
954 IncEvCountTotalESDs();
958 //____________________________________________________________________________
959 void AliTOFQADataMakerRec::StartOfDetectorCycle()
962 //Detector specific actions at start of cycle
963 fCalibData = GetCalibData();
968 //____________________________________________________________________________
969 void AliTOFQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
971 //Detector specific actions at end of cycle
972 // do the QA checking
973 ResetEventTrigClasses();
975 for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
976 if ( !AliQAv1::Instance()->IsEventSpecieSet(specie) ) continue ;
977 SetEventSpecie(AliRecoParam::ConvertIndex(specie));
979 for (int itc=-1;itc<GetNTrigClasses();itc++) { // RS: loop over eventual clones per trigger class
981 if (fEnableDqmShifterOpt) {
982 // RS: fetch the histograms for given trigger class
983 TObjArray& arrRW = *GetRawsDataOfTrigClass(itc);
985 // Help make the raw qa histogram easier to interpret for the DQM shifter
986 if (!arrRW[ 0] || !arrRW[ 5] || !arrRW[10] || !arrRW[15] || !arrRW[16] || !arrRW[17]) continue;
988 printf("=========>Processed %i physics raw of specie %s with TrigGlass %d\n",
989 GetEvCountCycleRaws(itc),AliRecoParam::GetEventSpecieName(specie), itc);
991 //Double_t monitorPeriodLength=fProcessedRawEventN*600*1E-9;//in s
994 //set minima and maxima to allow log scale
995 Double_t yTimeMax = ((TH1*)arrRW[5])->GetMaximum()*1.05;
996 Double_t yTotMax = ((TH1*)arrRW[10])->GetMaximum()*1.05;
997 fLineExpTimeMin->SetY2(yTimeMax);
998 fLineExpTimeMax->SetY2(yTimeMax);
999 fLineExpTotMin->SetY2(yTotMax);
1000 fLineExpTotMax->SetY2(yTotMax);
1002 for (Int_t j=0;j<18;j++){
1003 if ((j==0)||(j==5)||(j==10)||(j==15)||(j==16)||(j==17)) {
1004 TH1* htmp = (TH1*)arrRW[j];
1005 htmp->GetXaxis()->SetLabelOffset(0.005);
1006 htmp->GetXaxis()->SetLabelSize(0.05);
1007 htmp->GetXaxis()->SetTitleOffset(0.8);
1008 htmp->GetXaxis()->SetTitleSize(0.05);
1009 htmp->GetYaxis()->SetLabelOffset(0.005);
1010 htmp->GetYaxis()->SetLabelSize(0.06);
1011 htmp->GetYaxis()->SetTitleOffset(0.8);
1012 htmp->GetYaxis()->SetTitleSize(0.06);
1015 //make up for all histos
1016 for(Int_t j=0;j<5;j++) {
1017 TH1* htmp = (TH1*)arrRW[j];
1018 if (!htmp) continue;
1019 htmp->SetMarkerColor(kBlue);
1020 htmp->SetMarkerStyle(8);
1021 htmp->SetMarkerSize(0.7);
1023 for(Int_t j=5;j<15;j++) {
1024 TH1* htmp = (TH1*)arrRW[j];
1025 if (!htmp) continue;
1026 htmp->SetLineColor(kBlue);
1027 htmp->SetLineWidth(1);
1028 htmp->SetMarkerColor(kBlue);
1031 TH1* htmp = (TH1*)arrRW[15];
1032 htmp->SetLineColor(kBlue);
1033 htmp->SetLineWidth(1);
1034 htmp->SetMarkerStyle(8);
1035 htmp->SetMarkerSize(0.7);
1036 htmp->SetMarkerColor(kBlue);//Option("bar");
1038 TString title25 = Form("Map of hit pads according to CTTM numbering (Max Fired Macropad = %i)",fgCutNmaxFiredMacropad);
1039 TString title26 = Form("Map of hit macropads according to CTTM numbering (Max Fired Macropad = %i)",fgCutNmaxFiredMacropad);
1040 TString title27 = Form("#Deltat vs #Delta#Phi of hit macropads (Max Fired Macropad = %i)",fgCutNmaxFiredMacropad);
1042 if ( (htmp=(TH1*)arrRW[16]) ) htmp->SetOption("colz");
1043 if ( (htmp=(TH1*)arrRW[17]) ) htmp->SetOption("colz");
1044 if ( (htmp=(TH1*)arrRW[18]) ) htmp->SetOption("colz");
1045 if ( (htmp=(TH1*)arrRW[22]) ) htmp->SetOption("colz");
1046 if ( (htmp=(TH1*)arrRW[23]) ) htmp->SetOption("colz");
1047 if ( (htmp=(TH1*)arrRW[24]) ) htmp->SetOption("colz");
1048 if ( (htmp=(TH1*)arrRW[25]) ) {
1049 htmp->SetOption("colz");
1050 htmp->SetTitle(title25.Data());
1052 if ( (htmp=(TH1*)arrRW[26]) ) {
1053 htmp->SetOption("colz");
1054 htmp->SetTitle(title26.Data());
1056 if ( (htmp=(TH1*)arrRW[27]) ){
1057 htmp->SetOption("colz");
1058 htmp->SetTitle(title27.Data());
1061 }//END ENABLE DQM SHIFTER OPT
1062 } // RS: loop over trigger classes
1065 AliQAChecker::Instance()->Run(AliQAv1::kTOF, task, list) ;
1068 //____________________________________________________________________________
1069 void AliTOFQADataMakerRec::GetMapIndeces(const Int_t* const in , Int_t* out)
1072 //return appropriate indeces for the theta-phi map
1075 Int_t npadX = AliTOFGeometry::NpadX();
1076 Int_t npadZ = AliTOFGeometry::NpadZ();
1077 Int_t nStripA = AliTOFGeometry::NStripA();
1078 Int_t nStripB = AliTOFGeometry::NStripB();
1079 Int_t nStripC = AliTOFGeometry::NStripC();
1081 Int_t isector = in[0];
1082 Int_t iplate = in[1];
1083 Int_t istrip = in[2];
1084 Int_t ipadX = in[3];
1085 Int_t ipadZ = in[4];
1087 Int_t stripOffset = 0;
1093 stripOffset = nStripC;
1096 stripOffset = nStripC+nStripB;
1099 stripOffset = nStripC+nStripB+nStripA;
1102 stripOffset = nStripC+nStripB+nStripA+nStripB;
1105 AliDebug(1,Form("Wrong plate number in TOF (%d) !",iplate));
1108 Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
1109 Int_t phiindex=npadX*isector+ipadX+1;
1115 //---------------------------------------------------------------
1116 Int_t AliTOFQADataMakerRec::GetStripIndex(const Int_t * const in)
1118 /* return tof strip index between 0 and 91 */
1120 Int_t nStripA = AliTOFGeometry::NStripA();
1121 Int_t nStripB = AliTOFGeometry::NStripB();
1122 Int_t nStripC = AliTOFGeometry::NStripC();
1124 // Int_t isector = in[0];
1125 Int_t iplate = in[1];
1126 Int_t istrip = in[2];
1127 //Int_t ipadX = in[3];
1128 //Int_t ipadZ = in[4];
1130 Int_t stripOffset = 0;
1136 stripOffset = nStripC;
1139 stripOffset = nStripC+nStripB;
1142 stripOffset = nStripC+nStripB+nStripA;
1145 stripOffset = nStripC+nStripB+nStripA+nStripB;
1148 AliDebug(1,Form("Wrong plate number in TOF (%d) !",iplate));
1153 if (stripOffset<0 || stripOffset>92) return -1;
1155 return (stripOffset+istrip);
1158 //---------------------------------------------------------------
1159 Bool_t AliTOFQADataMakerRec::CheckVolumeID(const Int_t * const volumeID)
1162 //Checks volume ID validity
1164 for (Int_t j=0;j<5;j++){
1165 if (volumeID[j]<0) {
1166 AliDebug(1,Form("Invalid detector volume index for volumeID[%i]",j));
1173 //---------------------------------------------------------------
1174 Bool_t AliTOFQADataMakerRec::CheckEquipID(const Int_t * const equipmentID)
1177 //Checks equipment ID validity
1178 for (Int_t j=0;j<5;j++){
1179 if (equipmentID[j]<0) {
1180 AliDebug(1,Form("Invalid equipment volume index for equipmentID[%i]",j));
1186 //---------------------------------------------------------------
1187 Bool_t AliTOFQADataMakerRec::FilterLTMData(const Int_t * const equipmentID) const
1189 /*It returns kTRUE if data come from LTM.
1190 It thus filters trigger-related signals */
1192 Int_t ddl, trm, tdc;
1193 //if (!CheckEquipID(equipmentID)) return kFALSE;
1194 ddl = equipmentID[0];
1195 trm = equipmentID[1];
1196 tdc = equipmentID[3];
1198 if ((ddl%2==1) && (trm==3) && (tdc>11 && tdc<15))
1204 //---------------------------------------------------------------
1205 Bool_t AliTOFQADataMakerRec::FilterSpare(const Int_t * const equipmentID) const
1207 /*It returns kTRUE if data come from spare
1209 So far only check on TRM 3 crate left is implemented */
1211 Int_t ddl, trm, tdc;
1212 //if (!CheckEquipID(equipmentID)) return kFALSE;
1213 ddl = equipmentID[0];
1214 trm = equipmentID[1];
1215 tdc = equipmentID[3];
1217 if ((ddl%2==1) && (trm==3) && (tdc>2 && tdc<12))
1223 //-----------------------------------------------------------------------------
1224 void AliTOFQADataMakerRec::GetGeo2LTMIndex(const Int_t * const detind, Int_t *indexLTM) {
1226 // getting LTMmatrix indexes for current digit
1228 Int_t stripId=GetStripIndex(detind);
1230 if (detind[1]==0 || detind[1]==1 || (detind[1]==2 && detind[2]<=7)) { //A side
1231 if (detind[4]<24){ //R
1232 indexLTM[0] = detind[0]*2;
1234 indexLTM[0] = detind[0]*2+1;
1236 indexLTM[1]=stripId;
1240 indexLTM[0] = detind[0]*2+36;
1242 indexLTM[0] = (detind[0]*2+1)+36;
1244 indexLTM[1]=90-stripId;
1247 // if (indexLTM[0]<36) { //A side
1248 // if (detind[1] ==0){
1249 // indexLTM[1] = detind[2];
1251 // else if (detind[1] ==1){
1252 // indexLTM[1] = detind[2]+nStripB;
1254 // else if (detind[1] ==2){
1255 // indexLTM[1] = detind[2]+19*2;
1258 // AliError("Smth Wrong!!!");
1262 // if (detind[1]==2){
1263 // if (detind[4]<24)
1264 // indexLTM[1] = (nStripAL-detind[2])+nStripC+nStripB;
1266 // indexLTM[1] = (nStripAR-detind[2])+nStripC+nStripB;
1268 // else if (detind[1] ==3){
1269 // indexLTM[1] = (nStripB-detind[2])+nStripC;
1271 // else if (detind[1] ==4){
1272 // indexLTM[1] = nStripC-detind[2];
1275 // AliError("Smth Wrong!!!");
1280 //-----------------------------------------------------------------------------
1281 void AliTOFQADataMakerRec::GetGeo2CTTMIndex(Int_t *detind, Int_t *indexCTTM) {
1283 // Returns CTTM index corresponding to the detector element detind
1285 GetGeo2LTMIndex(detind,indexCTTM);
1290 //-------------------------------------------------------------------------
1291 Int_t AliTOFQADataMakerRec::GetNumberOfFiredMacropad(AliRawReader * rawReader){
1294 TClonesArray * clonesRawData;
1295 if (!rawReader) return 0;
1296 for (Int_t iDDL = 0; iDDL < AliTOFGeometry::NDDL()*AliTOFGeometry::NSectors(); iDDL++){
1298 fTOFRawStream.LoadRawDataBuffersV2(iDDL);
1299 clonesRawData = (TClonesArray*)fTOFRawStream.GetRawData();
1300 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
1301 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
1302 if (tofRawDatum->GetTOF()){
1303 if ( (tofRawDatum->GetTRM()==3)&&
1304 (tofRawDatum->GetTDC()>11)&&
1305 (tofRawDatum->GetTDC()<15)) {
1314 //----------------------------------------------------------------
1315 void AliTOFQADataMakerRec::GetCTTMIndex(Int_t *equipid, Int_t *indexCTTM) {
1317 // Returns CTTM index corresponding to the equipment id equipid, only for LTM hits
1318 // equipid = (crate, trm, chain, tdc, channel)
1320 if ( (equipid[1]!=3)||(equipid[3]<12) ){
1325 Int_t modStrip2LTM[3][8]={ { 0, 1, 2, 3, 4, 5, 6, 7},
1326 { 8, 9, 10, 11, 12, 13, 14, 15},
1327 {16, 17, 18, 19, 20, 21, 22, 23}
1330 Int_t DDL2LTMmatrix[72]={0,1,37,36,2,3,39,38,4,5,41,40,6,7,43,42,8,9,45,44,10,11,47,46,12,13,49,48,14,15,51,50,16,17,53,52,18,19,
1331 55,54,20,21,57,56,22,23,59,58,24,25,61,60,26,27,63,62,28,29,65,64,30,31,67,66,32,33,69,68,34,35,71,70};
1333 Int_t itdc=equipid[3]%12;
1337 else crate=equipid[0];
1339 indexCTTM[0]=DDL2LTMmatrix[crate];
1340 indexCTTM[1]=modStrip2LTM[itdc][equipid[4]];
1344 //_____________________________________________________________________________
1345 void AliTOFQADataMakerRec::ReadHistogramRangeFromFile(const Char_t * filename)
1348 // read histogram ranges from configuration file
1351 AliInfo("Config file with histograms ranges not found or invalid -> use default values.");
1352 SetDefaultHistogramRange();
1353 SetDefaultCutNmaxFiredMacropad();
1357 FILE * configFile = fopen(filename,"r");
1359 AliInfo("Cannot open config file with histograms ranges -> use default values.");
1360 SetDefaultHistogramRange();
1364 if (feof(configFile)){
1365 AliInfo("Unexpected EOF of config file with histograms ranges -> use default values.");
1366 SetDefaultHistogramRange();
1370 Int_t minMulti=9999, maxMulti=-9999;
1371 Int_t nbinsMulti=0,nbinsTime=0;
1372 Float_t minTime=9999.0, maxTime=-9999.0;
1373 Int_t cutFiredMacropad=0;
1375 fscanf(configFile,"%10i %10i %10i %10f %10f", &cutFiredMacropad,&minMulti,&maxMulti,&minTime,&maxTime);
1377 //set multiplicity histo ranges
1378 if (minMulti>maxMulti){
1379 AliInfo("Invalid range for multiplicity histogram set. Changing to defualt values.");
1380 SetDefaultMultiHistogramRange();
1382 nbinsMulti = maxMulti-minMulti;
1383 SetNbinsMultiplicityHisto(nbinsMulti);
1384 SetMultiplicityHistoRange(minMulti,maxMulti);
1385 AliInfo(Form("Setting multiplicity histogram ranges to: multMin = %i - multMax = %i - nMultBins = %i",
1386 fgRangeMinMultiplicity, fgRangeMaxMultiplicity, fgNbinsMultiplicity));
1389 //set time histo ranges
1390 if (minTime>maxTime){
1391 AliInfo("Invalid range for time histogram set. Changing to defualt values.");
1392 SetDefaultTimeHistogramRange();
1394 nbinsTime = TMath::Nint((maxTime - minTime)/fgkNbinsWidthTime);//ns
1395 maxTime=minTime+nbinsTime*fgkNbinsWidthTime;//ns
1396 SetNbinsTimeHisto(nbinsTime);
1397 SetTimeHistoRange(minTime,maxTime);
1398 AliInfo(Form("Setting time histogram ranges to: timeMin = %5.2f ns - timeMax = %5.2f ns - nTimeBins = %i",
1399 fgRangeMinTime, fgRangeMaxTime,fgNbinsTime));
1402 if ((cutFiredMacropad>0)&&(cutFiredMacropad<fgkFiredMacropadLimit)){
1403 AliInfo("Invalid value for cut on fired macropad. Changing to default values.");
1404 SetDefaultCutNmaxFiredMacropad();
1406 SetCutNmaxFiredMacropad(cutFiredMacropad);
1407 AliInfo(Form("Setting cut on fired macropad to: = %i",cutFiredMacropad));
1414 //_____________________________________________________________________________
1415 void AliTOFQADataMakerRec::SetDefaultHistogramRange()
1418 // set default histogram ranges (tuned on 2011 pp collisions)
1420 AliInfo("Setting all histogram ranges to default values.");
1421 SetDefaultMultiHistogramRange();
1422 SetDefaultTimeHistogramRange();
1423 SetDefaultCutNmaxFiredMacropad();
1427 //_____________________________________________________________________________
1428 void AliTOFQADataMakerRec::SetDefaultMultiHistogramRange()
1431 // set default histogram ranges (tuned on 2011 pp collisions)
1433 SetMultiplicityHistoRange (0, 200);
1434 SetNbinsMultiplicityHisto(200);
1435 AliInfo("Setting Multiplicity histogram ranges to default values.");
1436 AliInfo(Form("multMin = %i - multMax = %i - nMultBins = %i",
1437 fgRangeMinMultiplicity, fgRangeMaxMultiplicity, fgNbinsMultiplicity));
1441 //_____________________________________________________________________________
1442 void AliTOFQADataMakerRec::SetDefaultTimeHistogramRange()
1445 // set default histogram ranges (tuned on 2011 pp collisions)
1447 SetNbinsTimeHisto(250);
1448 SetTimeHistoRange (0.0,610.);
1450 AliInfo("Setting Time histogram ranges to default values:");
1451 AliInfo(Form("timeMin = %5.2f ns - timeMax = %5.2f ns - nTimeBins = %i",
1452 fgRangeMinTime, fgRangeMaxTime,fgNbinsTime));
1456 //------------------------------------------------------------------
1457 void AliTOFQADataMakerRec::SetDefaultCutNmaxFiredMacropad()
1460 // set default cut on fired macropad
1462 SetCutNmaxFiredMacropad(5);
1463 AliInfo(Form("Setting cut on fired macropad to default values: NfiredMacropad = %i", fgCutNmaxFiredMacropad));