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 18/01/2011
26 - reduced histo binning to reduce size
27 - added decoding errors plot
28 - added channel maps and options for DQM shifters
29 - new list of recPoints and ESDs plots
30 - removed hTrmChannels035 and hTrmChannels3671
32 Modified by bvonhall on 03/11/2010
33 - modified declaration of hTrmChannels035 and hTrmChannels3671 in EndOfDetectorCycle()
34 to prevent memory corruption
36 Modified by adecaro on 18/10/2010
37 - fTOFRawStream object set as private member
39 Modified by fbellini on 13/09/2010
40 - Set TLines as private members
41 - Set image flag for expert histos
43 Modified by fbellini on 14/06/2010
45 - use LoadRawDataBuffersV2()
47 Modified by fbellini on 10/05/2010
48 - Fixed EndOfDetectorCycle() memory corruption bug
50 Modified by fbellini on 22/04/2010
51 - Added filter for physics events
53 Modified by fbellini on 16/04/2010
54 - Added EnableDqmShifterOpt()
55 - Modified EndOfDetectorCycle() with options for DQM
58 Modified by fbellini on 30/03/2010
59 - Changed raws time histos range to 610ns
60 - Added FilterLTMData() and FilterSpare() methods
61 - Added check on enabled channels for raw data
62 - Updated RecPoints QA
64 Modified by fbellini on 02/03/2010
65 - Fixed raw data decoding methods (use AliTOFRawStream::LoadRawDataBuffer())
66 - Added filter for noisy channels and read map from OCDB
67 - Added GetCalibData() method
68 - Added CheckVolumeID() and CheckEquipID() methods
72 #include <TClonesArray.h>
76 #include <TPaveText.h>
78 #include "AliCDBManager.h"
79 #include "AliCDBEntry.h"
80 #include "AliESDEvent.h"
81 #include "AliESDtrack.h"
82 #include "AliQAChecker.h"
83 #include "AliRawReader.h"
84 #include "AliTOFRawStream.h"
85 #include "AliTOFcluster.h"
86 #include "AliTOFQADataMakerRec.h"
87 #include "AliTOFrawData.h"
88 #include "AliTOFGeometry.h"
89 #include "AliTOFChannelOnlineStatusArray.h"
90 #include "AliTOFDecoderSummaryData.h"
91 #include "AliTOFDecoderV2.h"
93 ClassImp(AliTOFQADataMakerRec)
95 //____________________________________________________________________________
96 AliTOFQADataMakerRec::AliTOFQADataMakerRec() :
97 AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kTOF), "TOF Quality Assurance Data Maker"),
99 fEnableNoiseFiltering(kFALSE),
100 fEnableDqmShifterOpt(kFALSE),
101 fProcessedRawEventN(0),
103 fLineExpTimeMin(new TLine(200., 0., 200., 0.)),
104 fLineExpTimeMax(new TLine(250., 0., 250., 0.)),
105 fLineExpTotMin(new TLine(5., 0., 5., 0.)),
106 fLineExpTotMax(new TLine(20., 0., 20., 0.)),
107 fTOFRawStream(AliTOFRawStream()),
108 fDecoderSummary(new AliTOFDecoderSummaryData())
113 for (Int_t sm=0;sm<17;sm++){
114 fLineSMid[sm] = new TLine( sm+1, 0., sm+1, 91.);
116 //initialize all TRM counters to -1 i.e. invalid value
117 // for (Int_t trm=0;trm<720;trm++){
118 // fTRMNoisyArray[trm]=-1;
119 // fTRMHwOkArray[trm]=-1;
120 // fTRMEnabledArray[trm]=-1;
125 //____________________________________________________________________________
126 AliTOFQADataMakerRec::AliTOFQADataMakerRec(const AliTOFQADataMakerRec& qadm) :
128 fCalibData(qadm.fCalibData),
129 fEnableNoiseFiltering(qadm.fEnableNoiseFiltering),
130 fEnableDqmShifterOpt(qadm.fEnableDqmShifterOpt),
131 fProcessedRawEventN(qadm.fProcessedRawEventN),
133 fLineExpTimeMin(qadm.fLineExpTimeMin),
134 fLineExpTimeMax(qadm.fLineExpTimeMax),
135 fLineExpTotMin(qadm.fLineExpTotMin),
136 fLineExpTotMax(qadm.fLineExpTotMax),
137 fTOFRawStream(qadm.fTOFRawStream),
138 fDecoderSummary(qadm.fDecoderSummary)
143 SetName((const char*)qadm.GetName()) ;
144 SetTitle((const char*)qadm.GetTitle());
146 for (Int_t sm=0;sm<17;sm++){
147 fLineSMid[sm]=qadm.fLineSMid[sm];
150 // for (Int_t trm=0;trm<10;trm++){
151 // fTRMNoisyArray[trm]=qadm.fTRMNoisyArray[trm];
152 // fTRMHwOkArray[trm]=qadm.fTRMHwOkArray[trm];
153 // fTRMEnabledArray[trm]=qadm.fTRMEnabledArray[trm];
158 //__________________________________________________________________
159 AliTOFQADataMakerRec& AliTOFQADataMakerRec::operator = (const AliTOFQADataMakerRec& qadm )
162 // assignment operator.
164 this->~AliTOFQADataMakerRec();
165 new(this) AliTOFQADataMakerRec(qadm);
169 //----------------------------------------------------------------------------
170 AliTOFQADataMakerRec::~AliTOFQADataMakerRec()
173 fTOFRawStream.Clear();
174 delete fLineExpTimeMin;
175 delete fLineExpTimeMax;
176 delete fLineExpTotMin;
177 delete fLineExpTotMax;
178 for (Int_t sm=0;sm<10;sm++){
179 delete fLineSMid[sm];
182 delete fDecoderSummary;
184 //----------------------------------------------------------------------------
185 AliTOFChannelOnlineStatusArray* AliTOFQADataMakerRec::GetCalibData() const
188 // Retrive TOF calib objects from OCDB
190 AliCDBManager *man = AliCDBManager::Instance();
193 cdbe = man->Get("TOF/Calib/Status",fRun);
195 AliWarning("Load of calibration data from default storage failed!");
196 AliWarning("Calibration data will be loaded from local storage ($ALICE_ROOT)");
197 man->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
198 cdbe = man->Get("TOF/Calib/Status",fRun);
200 // Retrieval of data in directory TOF/Calib/Data:
202 AliTOFChannelOnlineStatusArray * array = 0;
203 if (cdbe) array = (AliTOFChannelOnlineStatusArray *)cdbe->GetObject();
204 if (!array) AliFatal("No calibration data from calibration database !");
209 //____________________________________________________________________________
210 void AliTOFQADataMakerRec::InitRaws()
213 // create Raws histograms in Raws subdir
215 const Bool_t expert = kTRUE ;
216 const Bool_t saveCorr = kTRUE ;
217 const Bool_t image = kTRUE ;
219 TH1I * h0 = new TH1I("hTOFRaws", "TOF raw hit multiplicity; TOF raw hits number; Events ",200, 0, 200) ;
221 TH1I * h1 = new TH1I("hTOFRawsIA", "TOF raw hit multiplicity - I/A side; TOF raw hits number ;Events ",200, 0, 200) ;
222 TH1I * h2 = new TH1I("hTOFRawsOA", "TOF raw hit multiplicity - O/A side; TOF raw hits number ;Events ",200, 0, 200) ;
223 TH1I * h3 = new TH1I("hTOFRawsIC", "TOF raw hit multiplicity - I/C side; TOF raw hits number ;Events ",200, 0, 200) ;
224 TH1I * h4 = new TH1I("hTOFRawsOC", "TOF raw hit multiplicity - O/C side; TOF raw hits number ;Events ",200, 0, 200) ;
226 TH1F * h5 = new TH1F("hTOFRawsTime", "TOF Raws - Hit time (ns);Measured Hit time [ns];Hits", 250,0. ,610.) ;
227 TH1F * h6 = new TH1F("hTOFRawsTimeIA", "TOF Raws - Hit time (ns) - I/A side;Measured Hit time [ns];Hits", 250,0. ,610.) ;
228 TH1F * h7 = new TH1F("hTOFRawsTimeOA", "TOF Raws - Hit time (ns) - O/A side;Measured Hit time [ns];Hits", 250,0. ,610.) ;
229 TH1F * h8 = new TH1F("hTOFRawsTimeIC", "TOF Raws - Hit time (ns) - I/C side;Measured Hit time [ns];Hits", 250,0. ,610.) ;
230 TH1F * h9 = new TH1F("hTOFRawsTimeOC", "TOF Raws - Hit time (ns) - O/C side;Measured Hit time [ns];Hits", 250,0. ,610.) ;
232 TH1F * h10 = new TH1F("hTOFRawsToT", "TOF Raws - Hit ToT (ns);Measured Hit ToT (ns);Hits", 100, 0., 48.8) ;
234 TH1F * h11 = new TH1F("hTOFRawsToTIA", "TOF Raws - Hit ToT (ns) - I/A side;Measured Hit ToT (ns);Hits", 100, 0., 48.8) ;
235 TH1F * h12 = new TH1F("hTOFRawsToTOA", "TOF Raws - Hit ToT (ns) - O/A side;Measured Hit ToT (ns);Hits", 100, 0., 48.8) ;
236 TH1F * h13 = new TH1F("hTOFRawsToTIC", "TOF Raws - Hit ToT (ns) - I/C side;Measured Hit ToT (ns);Hits", 100, 0., 48.8) ;
237 TH1F * h14 = new TH1F("hTOFRawsToTOC", "TOF Raws - Hit ToT (ns) - O/C side;Measured Hit ToT (ns);Hits", 100, 0., 48.8) ;
239 TH1F * h15 = new TH1F("hTOFRawsLTMHits", "LTMs OR signals; Crate; Counts", 72, 0., 72.);
240 TH2F * h16 = new TH2F("hTOFrefMap", "TOF enabled channel reference map;sector;strip", 72, 0., 18., 91, 0., 91.);
241 TH2F * h17 = new TH2F("hTOFRawHitMap","TOF raw hit map;sector;strip", 72, 0., 18., 91, 0., 91.);
242 TH2I * h18 = new TH2I("hTOFDecodingErrors","Decoding error monitoring; DDL; Error ", 72, 0, 72, 13,1,14);
244 h18->GetYaxis()->SetBinLabel(1,"DRM ");
245 h18->GetYaxis()->SetBinLabel(2,"LTM ");
246 h18->GetYaxis()->SetBinLabel(3,"TRM 3 ");
247 h18->GetYaxis()->SetBinLabel(4,"TRM 4");
248 h18->GetYaxis()->SetBinLabel(5,"TRM 5");
249 h18->GetYaxis()->SetBinLabel(6,"TRM 6");
250 h18->GetYaxis()->SetBinLabel(7,"TRM 7");
251 h18->GetYaxis()->SetBinLabel(8,"TRM 8");
252 h18->GetYaxis()->SetBinLabel(9,"TRM 9");
253 h18->GetYaxis()->SetBinLabel(10,"TRM 10");
254 h18->GetYaxis()->SetBinLabel(11,"TRM 11");
255 h18->GetYaxis()->SetBinLabel(12,"TRM 12");
256 h18->GetYaxis()->SetBinLabel(13,"recovered");
258 TH1F * h19 = new TH1F("hTOFOrphansTime", "TOF Raws - Orphans time (ns);Measured Hit time [ns];Hits", 250, 0. ,610.) ;
259 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., 250, 0., 610.0) ;
260 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., 250, 0., 610.0) ;
261 TH2F * h22 = new TH2F("hTOFTimeVsStrip","TOF Raws - Hit time vs. strip (theta); Strip index;Raws TOF time (ns) ", 91,0.,91, 250, 0., 610.) ;
286 //add lines for DQM shifter
287 //fLineExpTimeMin = new TLine(200., 0., 200., 0.);
288 fLineExpTimeMin->SetLineColor(kGreen);
289 fLineExpTimeMin->SetLineWidth(2);
291 //fLineExpTimeMax = new TLine(250., 0., 250., 0.);
292 fLineExpTimeMax->SetLineColor(kGreen);
293 fLineExpTimeMax->SetLineWidth(2);
295 //fLineExpTotMin = new TLine( 5., 0., 5., 0.);
296 fLineExpTotMin->SetLineColor(kGreen);
297 fLineExpTotMin->SetLineWidth(2);
299 //fLineExpTotMax = new TLine(20., 0., 20., 0.);
300 fLineExpTotMax->SetLineColor(kGreen);
301 fLineExpTotMax->SetLineWidth(2);
303 for (Int_t sm=0;sm<17;sm++){
304 //fLineSMid[sm] = new TLine( 1+sm, 0., 1+sm, 91.);
305 fLineSMid[sm]->SetLineColor(kMagenta);
306 fLineSMid[sm]->SetLineWidth(2);
309 h5->GetListOfFunctions()->Add(fLineExpTimeMin);
310 h5->GetListOfFunctions()->Add(fLineExpTimeMax);
311 h10->GetListOfFunctions()->Add(fLineExpTotMin);
312 h10->GetListOfFunctions()->Add(fLineExpTotMax);
314 for (Int_t sm=0;sm<17;sm++){
315 h16->GetListOfFunctions()->Add(fLineSMid[sm]);
316 h17->GetListOfFunctions()->Add(fLineSMid[sm]);
320 TPaveText *phosHoleBox=new TPaveText(13,38,16,53,"b");
321 phosHoleBox->SetFillStyle(0);
322 phosHoleBox->SetFillColor(kWhite);
323 phosHoleBox->SetLineColor(kMagenta);
324 phosHoleBox->SetLineWidth(2);
325 phosHoleBox->AddText("PHOS");
326 h16->GetListOfFunctions()->Add(phosHoleBox);
327 h17->GetListOfFunctions()->Add(phosHoleBox);
329 // h0->SetDrawOption("logy");
330 // h5->SetDrawOption("logy");
331 // h10->SetDrawOption("logy");
333 Add2RawsList(h0, 0, !expert, image, !saveCorr) ;
334 Add2RawsList(h1, 1, expert, !image, !saveCorr) ;
335 Add2RawsList(h2, 2, expert, !image, !saveCorr) ;
336 Add2RawsList(h3, 3, expert, !image, !saveCorr) ;
337 Add2RawsList(h4, 4, expert, !image, !saveCorr) ;
338 Add2RawsList(h5, 5, !expert, image, !saveCorr) ;
339 Add2RawsList(h6, 6, expert, !image, !saveCorr) ;
340 Add2RawsList(h7, 7, expert, !image, !saveCorr) ;
341 Add2RawsList(h8, 8, expert, !image, !saveCorr) ;
342 Add2RawsList(h9, 9, expert, !image, !saveCorr) ;
343 Add2RawsList(h10, 10, !expert, image, !saveCorr) ;
344 Add2RawsList(h11, 11, expert, !image, !saveCorr) ;
345 Add2RawsList(h12, 12, expert, !image, !saveCorr) ;
346 Add2RawsList(h13, 13, expert, !image, !saveCorr) ;
347 Add2RawsList(h14, 14, expert, !image, !saveCorr) ;
348 Add2RawsList(h15, 15, !expert, image, !saveCorr) ;
349 Add2RawsList(h16, 16, !expert, image, !saveCorr) ;
350 Add2RawsList(h17, 17, !expert, image, !saveCorr) ;
351 Add2RawsList(h18, 18, expert, !image, !saveCorr) ;
352 Add2RawsList(h19, 19, expert, !image, !saveCorr) ;
353 Add2RawsList(h20, 20, expert, !image, !saveCorr) ;
354 Add2RawsList(h21, 21, expert, !image, !saveCorr) ;
355 Add2RawsList(h22, 22, expert, !image, !saveCorr) ;
358 //____________________________________________________________________________
359 void AliTOFQADataMakerRec::InitRecPoints()
362 // create RecPoints histograms in RecPoints subdir
365 const Bool_t expert = kTRUE ;
366 const Bool_t image = kTRUE ;
368 TH1I * h0 = new TH1I("hTOFRecPoints", "TOF RecPoints multiplicity ; TOF RecPoints number;Events",200, 0, 200) ;
370 TH1F * h1 = new TH1F("hTOFRecPointsTimeIA", "RecPoints Time Spectrum in TOF (ns)- I/A side; Calibrated TOF time [ns];Events", 250, 0., 610.) ;
371 TH1F * h2 = new TH1F("hTOFRecPointsTimeOA", "RecPoints Time Spectrum in TOF (ns)- O/A side; Calibrated TOF time [ns];Events", 250, 0., 610.) ;
372 TH1F * h3 = new TH1F("hTOFRecPointsTimeIC", "RecPoints Time Spectrum in TOF (ns)- I/C side; Calibrated TOF time [ns];Events", 250, 0., 610.) ;
373 TH1F * h4 = new TH1F("hTOFRecPointsTimeOC", "RecPoints Time Spectrum in TOF (ns)- O/C side; Calibrated TOF time [ns];Events", 250, 0., 610.) ;
375 TH1F * h5 = new TH1F("hTOFRecPointsRawTimeIA", "RecPoints raw Time Spectrum in TOF (ns)-I/A side; Measured TOF time [ns];Hits", 250, 0., 610.) ;
376 TH1F * h6 = new TH1F("hTOFRecPointsRawTimeOA", "RecPoints raw Time Spectrum in TOF (ns)-O/A side; Measured TOF time [ns];Hits", 250, 0., 610.) ;
377 TH1F * h7 = new TH1F("hTOFRecPointsRawTimeIC", "RecPoints raw Time Spectrum in TOF (ns)-I/C side; Measured TOF time [ns];Hits", 250, 0., 610.) ;
378 TH1F * h8 = new TH1F("hTOFRecPointsRawTimeOC", "RecPoints raw Time Spectrum in TOF (ns)-O/C side; Measured TOF time [ns];Hits", 250, 0., 610.) ;
380 TH1F * h9 = new TH1F("hTOFRecPointsToTIA", "RecPoints ToT Spectrum in TOF (ns)-I/A side; Measured TOT [ns];Hits", 100, 0., 48.8 ) ;
381 TH1F * h10 = new TH1F("hTOFRecPointsToTOA", "RecPoints ToT Spectrum in TOF (ns)-O/A side; Measured TOT [ns];Hits", 100, 0., 48.8 ) ;
382 TH1F * h11 = new TH1F("hTOFRecPointsToTIC", "RecPoints ToT Spectrum in TOF (ns)-I/C side; Measured TOT [ns];Hits", 100, 0., 48.8 ) ;
383 TH1F * h12 = new TH1F("hTOFRecPointsToTOC", "RecPoints ToT Spectrum in TOF (ns)-O/C side; Measured TOT [ns];Hits", 100, 0., 48.8 ) ;
385 TH2F * h13 = new TH2F("hTOFRecPointsClusMap","RecPoints map; sector ;strip", 72, 0., 18., 91, 0., 91.) ;
386 TH2F * h14 = new TH2F("hTOFRecPointsTimeVsStrip","RecPoints TOF time vs. strip (theta); Strip index; RecPoints TOF time (ns) ",91, 0., 91., 250, 0., 610.) ;
387 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) ;
388 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) ;
408 Add2RecPointsList(h0, 0, !expert, image) ;
409 Add2RecPointsList(h1, 1, !expert, image) ;
410 Add2RecPointsList(h2, 2, !expert, image) ;
411 Add2RecPointsList(h3, 3, !expert, image) ;
412 Add2RecPointsList(h4, 4, !expert, image) ;
413 Add2RecPointsList(h5, 5, expert, !image) ;
414 Add2RecPointsList(h6, 6, expert, !image) ;
415 Add2RecPointsList(h7, 7, expert, !image) ;
416 Add2RecPointsList(h8, 8, expert, !image) ;
417 Add2RecPointsList(h9, 9, !expert, !image) ;
418 Add2RecPointsList(h10, 10, !expert, !image) ;
419 Add2RecPointsList(h11, 11, !expert, !image) ;
420 Add2RecPointsList(h12, 12, !expert, !image) ;
421 Add2RecPointsList(h13, 13, expert, !image) ;
422 Add2RecPointsList(h14, 14, expert, image) ;
423 Add2RecPointsList(h15, 15, expert, !image) ;
424 Add2RecPointsList(h16, 16, expert, !image) ;
428 //____________________________________________________________________________
429 void AliTOFQADataMakerRec::InitESDs()
432 //create ESDs histograms in ESDs subdir
435 const Bool_t expert = kTRUE ;
436 const Bool_t image = kTRUE ;
438 TH1I * h0 = new TH1I("hTOFESDs", "Number of matched TOF tracks per event;Number of TOF matched ESD tracks;Counts", 200, 0, 200) ;
439 TH1F * h1 = new TH1F("hTOFESDsTime", "Matched ESDs tracks: TOF Time spectrum; Calibrated TOF time [ns];Counts", 250, 0., 610. ) ;
440 TH1F * h2 = new TH1F("hTOFESDsRawTime", "Matched ESDs tracks: TOF raw Time spectrum;Measured TOF time [ns];Counts", 250, 0., 610.) ;
441 TH1F * h3 = new TH1F("hTOFESDsToT", "Matched ESDs tracks: TOF ToT spectrum; ESDs ToT [ns];Counts",100, 0., 48.8) ;
442 TH1F * h4 = new TH1F("hTOFESDskTOFOUT", "p_{T} distribution of tracks with kTOFout; p_{T} (GeV/c);Counts", 50, 0.20, 5.00) ;
443 TH1F * h5 = new TH1F("hTOFESDskTIME", "p_{T} distribution of tracks with kTOFout && kTIME; p_{T} (GeV/c);Counts", 50, 0.20, 5.00) ;
444 TH1F * h6 = new TH1F("hTOFESDsMatched", "p_{T} distribution of tracks with kTOFout && TOFtime>0; p_{T} (GeV/c);Counts", 50, 0.20, 5.00) ;
445 TH1F * h7 = new TH1F("hTOFESDsMatchingProb", "TPC-TOF track-matching probability;TOF matching probability (%) ;Counts",101, -1.0, 100) ;
446 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.) ;
447 TH1F * h9 = new TH1F("hTOFHitsLength", "Matched ESDs tracks: Length Spectrum; Track length [cm];Counts", 800, 0., 800) ;
460 Add2ESDsList(h0, 0, !expert, image) ;
461 Add2ESDsList(h1, 1, !expert, image) ;
462 Add2ESDsList(h2, 2, expert, !image) ;
463 Add2ESDsList(h3, 3, !expert, !image) ;
464 Add2ESDsList(h4, 4, expert, image) ;
465 Add2ESDsList(h5, 5, expert, image) ;
466 Add2ESDsList(h6, 6, expert, image) ;
467 Add2ESDsList(h7, 7, expert, image) ;
468 Add2ESDsList(h8, 8, expert, !image) ;
469 Add2ESDsList(h9, 9, !expert, !image) ;
474 //____________________________________________________________________________
475 void AliTOFQADataMakerRec::MakeRaws(AliRawReader* rawReader)
478 // makes data from Raws
480 if (rawReader->GetType()==7) {
482 Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;//in ns
483 Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
484 Int_t ntof[5]; /* 0=tot, 1=IA, 2=OA, 3=IC, 4=OC*/
485 for (Int_t j=0;j<5;j++){ ntof[j]=0;}
486 Int_t equipmentID[5]; //(ddl, trm, chain,tdc,channel)
487 Int_t volumeID[5]; //(sector,plate,strip,padX,padZ)
488 Int_t volumeID2[5]; //(sector,plate,strip,padZ,padX) to use AliTOFGeometry::GetIndex()
491 TClonesArray * clonesRawData;
492 fTOFRawStream.SetRawReader(rawReader);
494 //uncomment if needed to apply DeltaBC correction
495 //fTOFRawStream.ApplyBCCorrections(kTRUE);
497 fDecoderSummary->Reset();
498 for (Int_t iDDL = 0; iDDL < AliTOFGeometry::NDDL()*AliTOFGeometry::NSectors(); iDDL++){
500 fTOFRawStream.LoadRawDataBuffersV2(iDDL);
502 //* get decoding error counters
503 fDecoderSummary = ( (AliTOFDecoderV2*) fTOFRawStream.GetDecoderV2() )->GetDecoderSummaryData();
504 if ( (fDecoderSummary) && (fDecoderSummary ->GetErrorDetected()) ) {
505 Int_t errorSlotID=(Int_t) fDecoderSummary->GetErrorSlotID();
506 GetRawsData(18)->Fill(iDDL,errorSlotID);
507 if (fDecoderSummary -> GetRecoverError() )
508 GetRawsData(18)->Fill(iDDL,13);
511 clonesRawData = (TClonesArray*)fTOFRawStream.GetRawData();
512 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
513 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
515 if (tofRawDatum->GetTOF()){
517 equipmentID[1]=tofRawDatum->GetTRM();
518 equipmentID[2]=tofRawDatum->GetTRMchain();
519 equipmentID[3]=tofRawDatum->GetTDC();
520 equipmentID[4]=tofRawDatum->GetTDCchannel();
522 if (CheckEquipID(equipmentID)){
523 fTOFRawStream.EquipmentId2VolumeId(iDDL,
524 tofRawDatum->GetTRM(),
525 tofRawDatum->GetTRMchain(),
526 tofRawDatum->GetTDC(),
527 tofRawDatum->GetTDCchannel(),
530 if (FilterLTMData(equipmentID)) { //counts LTM hits
531 if (equipmentID[2]==1) { //crate left, A-side or C-side
532 GetRawsData(15)->Fill(equipmentID[0]);
534 if (equipmentID[0]<36) { GetRawsData(15)->Fill(equipmentID[0]-1); }
535 else { GetRawsData(15)->Fill(equipmentID[0]+1); }
541 if (CheckVolumeID(volumeID)){
542 volumeID2[0]=volumeID[0];
543 volumeID2[1]=volumeID[1];
544 volumeID2[2]=volumeID[2];
545 volumeID2[3]=volumeID[4];
546 volumeID2[4]=volumeID[3];
547 chIndex=AliTOFGeometry::GetIndex(volumeID2);
549 if (tofRawDatum->GetTOT()){
550 if (!(fCalibData->GetNoiseStatus(chIndex)==AliTOFChannelOnlineStatusArray::kTOFNoiseBad)
551 && (fCalibData->GetHWStatus(chIndex) == AliTOFChannelOnlineStatusArray::kTOFHWOk)) {//noise and enabled filter
552 ntof[0]++; //counter for tof hits
554 //fill global spectra for DQM plots
555 GetRawsData(5)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;//in ns
556 GetRawsData(10)->Fill( tofRawDatum->GetTOT()*tot2ns) ;//in ns
558 //fill side-related spectra for experts plots
559 Int_t ddlACside=iDDL/36; // 0 or 1
560 Int_t ddlPerSm=iDDL%4;
562 if (volumeID2[0]>4 && volumeID2[0]<14){ //O side
563 if (ddlPerSm<2){ //A side
565 GetRawsData(6)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;
566 GetRawsData(11)->Fill( tofRawDatum->GetTOT()*tot2ns) ;
569 GetRawsData(8)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;
570 GetRawsData(13)->Fill( tofRawDatum->GetTOT()*tot2ns) ;
573 if (volumeID2[0]<5 || volumeID2[0]>13){ //I side
574 if (ddlPerSm<2){ //A side
576 GetRawsData(7)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;
577 GetRawsData(12)->Fill( tofRawDatum->GetTOT()*tot2ns) ;
580 GetRawsData(9)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;
581 GetRawsData(14)->Fill( tofRawDatum->GetTOT()*tot2ns) ;
587 Int_t trm= iDDL*10+(equipmentID[1]-3);
588 GetRawsData(20+ddlACside)->Fill(trm,tofRawDatum->GetTOF()*tdc2ns);
589 GetRawsData(22)->Fill(GetStripIndex(volumeID),tofRawDatum->GetTOF()*tdc2ns) ;
590 Short_t fea = volumeID2[4]/12;
591 Float_t hitmapx = volumeID2[0] + ((Double_t)(3 - fea) + 0.5) *0.25;
592 GetRawsData(17)->Fill(hitmapx,GetStripIndex(volumeID2));
596 if (!(fCalibData->GetNoiseStatus(chIndex) == AliTOFChannelOnlineStatusArray::kTOFNoiseBad)
597 && (fCalibData->GetHWStatus(chIndex) == AliTOFChannelOnlineStatusArray::kTOFHWOk))
598 GetRawsData(19)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;//in ns
600 }//end volumeID check
604 clonesRawData->Clear();
607 for (Int_t j=0;j<5;j++) { GetRawsData(j)->Fill(ntof[j]); }
608 fProcessedRawEventN++;
609 fTOFRawStream.Clear();
611 AliDebug(1,Form("Event of type %d found. Skipping non-physics event for QA.\n", rawReader->GetType()));
614 //fill reference map for DQM shifter only once in a detector cycle
616 Int_t geoId[5]={-1,-1,-1,-1,-1};// pgeoId=(sm, mod, strip, padZ, padX)
617 Int_t detId[5]={-1,-1,-1,-1,-1};//detID=(ddl,trm,tdc, chain,channel)
619 for (Int_t ch = 0; ch < fCalibData->GetSize(); ch++) {
620 AliTOFGeometry::GetVolumeIndices(ch,geoId);
621 AliTOFRawStream::Geant2EquipmentId(geoId,detId);
622 if ((detId[1]<0)||(detId[0]<0)) continue;
623 trmIndex=(detId[1]-3)+detId[0]*10;
626 // if (fCalibData->GetNoiseStatus(ch)==AliTOFChannelOnlineStatusArray::kTOFNoiseBad)
627 // fTRMNoisyArray[trmIndex]+=1;
628 // if (fCalibData->GetHWStatus(ch) == AliTOFChannelOnlineStatusArray::kTOFHWOk)
629 // fTRMHwOkArray[trmIndex]+=1;
631 if ( (!(fCalibData->GetNoiseStatus(ch)==AliTOFChannelOnlineStatusArray::kTOFNoiseBad))
632 && (fCalibData->GetHWStatus(ch) == AliTOFChannelOnlineStatusArray::kTOFHWOk) ){
633 //fTRMEnabledArray[trmIndex]+=1;
634 //fill reference map with info from OCDB
635 Short_t fea = geoId[4]/12;
636 Float_t hitmapx = geoId[0] + ((Double_t)(3 - fea) + 0.5)*0.25;
637 GetRawsData(16)->Fill(hitmapx, GetStripIndex(geoId));
640 //printf("Counters for noisy, enabled and good channels in TOF TRMs read from OCDB.\n");
644 //enable options for DQM shifter
645 EnableDqmShifterOpt(kTRUE);
648 //____________________________________________________________________________
649 void AliTOFQADataMakerRec::MakeRecPoints(TTree * clustersTree)
652 // Make data from Clusters
655 Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
656 Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
658 Int_t volumeID2[5];//(sm, plate,strip, padZ,padX)
659 //Int_t volumeID[5];//(sm, plate,strip, padX,padZ)
661 TBranch *branch=clustersTree->GetBranch("TOF");
663 AliError("can't get the branch with the TOF clusters !");
667 static TClonesArray dummy("AliTOFcluster",10000);
669 TClonesArray *clusters=&dummy;
670 branch->SetAddress(&clusters);
673 clustersTree->GetEvent(0);
675 GetRecPointsData(0)->Fill((Int_t)clusters->GetEntriesFast()) ;
677 TIter next(clusters) ;
679 while ( (c = dynamic_cast<AliTOFcluster *>(next())) ) {
681 // volumeID2[0] = c->GetDetInd(0);
682 // volumeID2[1] = c->GetDetInd(1);
683 // volumeID2[2] = c->GetDetInd(2);
684 // volumeID2[3] = c->GetDetInd(4); //padX
685 // volumeID2[4] = c->GetDetInd(3); //padZ
687 for (Int_t i=0;i<5;i++){
688 volumeID2[i]=c->GetDetInd(i); //X and Z indeces inverted in RecPoints
690 //Int_t chIndex=AliTOFGeometry::GetIndex(volumeID2);
691 Int_t iDDL=AliTOFRawStream::Geant2DDL(volumeID2);
692 Int_t iTRM=AliTOFRawStream::Geant2TRM(volumeID2);
693 Short_t fea = volumeID2[4]/12;
694 Float_t hitmapx = volumeID2[0] + ((Double_t)(3 - fea) + 0.5) *0.25;
695 Int_t ddlACside=iDDL/36; // 0 or 1
696 Int_t ddlPerSm=iDDL%4;
698 if ((c->GetTDCRAW()) && (c->GetTDC()) && (c->GetToT())){
699 if (volumeID2[0]>4 && volumeID2[0]<14){ //I side
700 if (ddlPerSm<2){ //A side
701 GetRecPointsData(1)->Fill( c->GetTDC()*tdc2ns) ;//in ns
702 GetRecPointsData(5)->Fill( c->GetTDCRAW()*tdc2ns) ;//in ns
703 GetRecPointsData(9)->Fill( c->GetToT()*tot2ns) ;//in ns
705 GetRecPointsData(3)->Fill( c->GetTDC()*tdc2ns) ;//in ns
706 GetRecPointsData(7)->Fill( c->GetTDCRAW()*tdc2ns) ;//in ns
707 GetRecPointsData(11)->Fill( c->GetToT()*tot2ns) ;//in ns
710 if (volumeID2[0]<5 || volumeID2[0]>13){ //O side
711 if (ddlPerSm<2){ //A side
712 GetRecPointsData(2)->Fill( c->GetTDC()*tdc2ns) ;//in ns
713 GetRecPointsData(6)->Fill( c->GetTDCRAW()*tdc2ns) ;//in ns
714 GetRecPointsData(10)->Fill( c->GetToT()*tot2ns) ;//in ns
716 GetRecPointsData(4)->Fill( c->GetTDC()*tdc2ns) ;//in ns
717 GetRecPointsData(8)->Fill( c->GetTDCRAW()*tdc2ns) ;//in ns
718 GetRecPointsData(12)->Fill( c->GetToT()*tot2ns) ;//in ns
722 GetRecPointsData(13)->Fill(hitmapx,GetStripIndex(volumeID2));
723 GetRecPointsData(14)->Fill(GetStripIndex(volumeID2), c->GetTDC()*tdc2ns) ;
724 Int_t trm= iDDL*10+(iTRM-3);
725 if (ddlACside==0) { //A side
726 GetRecPointsData(15)->Fill(trm,c->GetTDC()*tdc2ns);
728 GetRecPointsData(16)->Fill(trm,c->GetTDC()*tdc2ns);
732 EnableDqmShifterOpt(kFALSE);
735 //____________________________________________________________________________
736 void AliTOFQADataMakerRec::MakeESDs(AliESDEvent * esd)
739 // make QA data from ESDs
741 const Double_t speedOfLight = TMath::C()*1E2*1E-12; // cm/ps
742 const Double_t pionMass = 0.13957018; //GeV/c^2
744 Int_t ntrk = esd->GetNumberOfTracks() ;
749 AliESDtrack *track=esd->GetTrack(ntrk);
750 Double_t tofTime=track->GetTOFsignal();//in ps
751 Double_t tofTimeRaw=track->GetTOFsignalRaw();//in ps
752 Double_t tofToT=track->GetTOFsignalToT(); //in ps
754 UInt_t status=track->GetStatus();
755 if (track->IsOn(AliESDtrack::kTPCrefit)) {
757 Double_t y=track->Eta();
758 if (TMath::Abs(y)<0.9) { //select TOF acceptance
759 if ((status&AliESDtrack::kTOFout)!=0) { //define matching
761 GetESDsData(1)->Fill(tofTime*1E-3);
762 GetESDsData(2)->Fill(tofTimeRaw*1E-3);
763 GetESDsData(3)->Fill(tofToT*1E-3);
764 GetESDsData(4)->Fill(track->Pt());
766 Double_t length =track->GetIntegratedLength();
767 Double_t mom2=(track->Pt()*track->Pt())+(track->Pz()*track->Pz());
768 Double_t piTexp = TMath::Sqrt(1+(pionMass*pionMass/mom2))*length/speedOfLight; //in ps
769 GetESDsData(8)->Fill(tofTime-piTexp);
770 GetESDsData(9)->Fill(length);
772 if ((status&AliESDtrack::kTIME)!=0)
773 GetESDsData(5)->Fill(track->Pt());
776 GetESDsData(6)->Fill(track->Pt());
777 } //end check on matched tracks
779 }//end check on TPCrefit
782 GetESDsData(0)->Fill(ntofout) ;
784 Float_t ratio = (Float_t)ntofout/(Float_t)ntpc*100.; //matching probability
785 GetESDsData(7)->Fill(ratio) ;
789 Float_t ratio = (Float_t)ntofout/(Float_t)ntpc*100; //matched over propagated to TOF outer radius
790 GetESDsData(8)->Fill(ratio) ;
792 EnableDqmShifterOpt(kFALSE);
795 //____________________________________________________________________________
796 void AliTOFQADataMakerRec::StartOfDetectorCycle()
799 //Detector specific actions at start of cycle
800 // ResetAllTRMcounters();
801 fCalibData = GetCalibData();
806 //____________________________________________________________________________
807 void AliTOFQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
809 //Detector specific actions at end of cycle
810 // do the QA checking
812 for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
813 if ( !AliQAv1::Instance()->IsEventSpecieSet(specie) )
816 if (fEnableDqmShifterOpt){
817 // Help make the raw qa histogram easier to interpret for the DQM shifter
818 if (!GetRawsData(0) || !GetRawsData(5) || !GetRawsData(10)
819 || !GetRawsData(15) || !GetRawsData(16) || !GetRawsData(17)) {
820 printf("No histogram for DQM found - Possible memory corruption ???. Please check\n") ;
823 printf("=========>Processed %i physics raw events. \n",fProcessedRawEventN);
825 //Double_t monitorPeriodLength=fProcessedRawEventN*600*1E-9;//in s
828 //set minima and maxima to allow log scale
829 Double_t yTimeMax = GetRawsData(5)->GetMaximum()*1.05;
830 Double_t yTotMax = GetRawsData(10)->GetMaximum()*1.05;
831 fLineExpTimeMin->SetY2(yTimeMax);
832 fLineExpTimeMax->SetY2(yTimeMax);
833 fLineExpTotMin->SetY2(yTotMax);
834 fLineExpTotMax->SetY2(yTotMax);
836 for (Int_t j=0;j<18;j++){
837 if ((j==0)||(j==5)||(j==10)||(j==15)||(j==16)||(j==17)) {
838 GetRawsData(j)->GetXaxis()->SetLabelOffset(0.005);
839 GetRawsData(j)->GetXaxis()->SetLabelSize(0.05);
840 GetRawsData(j)->GetXaxis()->SetTitleOffset(0.8);
841 GetRawsData(j)->GetXaxis()->SetTitleSize(0.05);
842 GetRawsData(j)->GetYaxis()->SetLabelOffset(0.005);
843 GetRawsData(j)->GetYaxis()->SetLabelSize(0.06);
844 GetRawsData(j)->GetYaxis()->SetTitleOffset(0.8);
845 GetRawsData(j)->GetYaxis()->SetTitleSize(0.06);
848 //make up for all histos
849 for(Int_t j=0;j<5;j++){
850 GetRawsData(j)->SetMarkerColor(kBlue);
851 GetRawsData(j)->SetMarkerStyle(8);
852 GetRawsData(j)->SetMarkerSize(0.7);
854 for(Int_t j=5;j<15;j++){
855 GetRawsData(j)->SetLineColor(kBlue);
856 GetRawsData(j)->SetLineWidth(1);
857 GetRawsData(j)->SetMarkerColor(kBlue);
858 //GetRawsData(j)->SetFillColor(kWhite);
859 //GetRawsData(j)->SetDrawOption("bar");
862 GetRawsData(15)->SetLineColor(kBlue);
863 GetRawsData(15)->SetLineWidth(1);
864 GetRawsData(15)->SetMarkerStyle(8);
865 GetRawsData(15)->SetMarkerSize(0.7);
866 GetRawsData(15)->SetMarkerColor(kBlue);//Option("bar");
868 GetRawsData(16)->SetOption("colz");
869 GetRawsData(17)->SetOption("colz");
870 GetRawsData(18)->SetOption("colz");
872 }//END ENABLE DQM SHIFTER OPT
874 AliQAChecker::Instance()->Run(AliQAv1::kTOF, task, list) ;
876 //____________________________________________________________________________
877 void AliTOFQADataMakerRec::GetMapIndeces(const Int_t* const in , Int_t* out)
880 //return appropriate indeces for the theta-phi map
883 Int_t npadX = AliTOFGeometry::NpadX();
884 Int_t npadZ = AliTOFGeometry::NpadZ();
885 Int_t nStripA = AliTOFGeometry::NStripA();
886 Int_t nStripB = AliTOFGeometry::NStripB();
887 Int_t nStripC = AliTOFGeometry::NStripC();
889 Int_t isector = in[0];
890 Int_t iplate = in[1];
891 Int_t istrip = in[2];
895 Int_t stripOffset = 0;
901 stripOffset = nStripC;
904 stripOffset = nStripC+nStripB;
907 stripOffset = nStripC+nStripB+nStripA;
910 stripOffset = nStripC+nStripB+nStripA+nStripB;
913 AliDebug(1,Form("Wrong plate number in TOF (%d) !",iplate));
916 Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
917 Int_t phiindex=npadX*isector+ipadX+1;
923 //---------------------------------------------------------------
924 Int_t AliTOFQADataMakerRec::GetStripIndex(const Int_t * const in)
926 /* return tof strip index between 0 and 91 */
928 Int_t nStripA = AliTOFGeometry::NStripA();
929 Int_t nStripB = AliTOFGeometry::NStripB();
930 Int_t nStripC = AliTOFGeometry::NStripC();
932 // Int_t isector = in[0];
933 Int_t iplate = in[1];
934 Int_t istrip = in[2];
935 //Int_t ipadX = in[3];
936 //Int_t ipadZ = in[4];
938 Int_t stripOffset = 0;
944 stripOffset = nStripC;
947 stripOffset = nStripC+nStripB;
950 stripOffset = nStripC+nStripB+nStripA;
953 stripOffset = nStripC+nStripB+nStripA+nStripB;
956 AliDebug(1,Form("Wrong plate number in TOF (%d) !",iplate));
961 if (stripOffset<0 || stripOffset>92) return -1;
963 return (stripOffset+istrip);
966 //---------------------------------------------------------------
967 Bool_t AliTOFQADataMakerRec::CheckVolumeID(const Int_t * const volumeID)
970 //Checks volume ID validity
973 for (Int_t j=0;j<5;j++){
975 AliDebug(1,Form("Invalid detector volume index for volumeID[%i]",j));
983 //---------------------------------------------------------------
984 Bool_t AliTOFQADataMakerRec::CheckEquipID(const Int_t * const equipmentID)
987 //Checks equipment ID validity
989 for (Int_t j=0;j<5;j++){
990 if (equipmentID[j]<0) {
991 AliDebug(1,Form("Invalid equipment volume index for equipmentID[%i]",j));
997 //---------------------------------------------------------------
998 Bool_t AliTOFQADataMakerRec::FilterLTMData(const Int_t * const equipmentID) const
1000 /*It returns kTRUE if data come from LTM.
1001 It thus filters trigger-related signals */
1003 Int_t ddl, trm, tdc;
1004 //if (!CheckEquipID(equipmentID)) return kFALSE;
1005 ddl = equipmentID[0];
1006 trm = equipmentID[1];
1007 tdc = equipmentID[3];
1009 if ((ddl%2==1) && (trm==3) && (tdc>11 && tdc<15))
1015 //---------------------------------------------------------------
1016 Bool_t AliTOFQADataMakerRec::FilterSpare(const Int_t * const equipmentID) const
1018 /*It returns kTRUE if data come from spare
1020 So far only check on TRM 3 crate left is implemented */
1022 Int_t ddl, trm, tdc;
1023 //if (!CheckEquipID(equipmentID)) return kFALSE;
1024 ddl = equipmentID[0];
1025 trm = equipmentID[1];
1026 tdc = equipmentID[3];
1028 if ((ddl%2==1) && (trm==3) && (tdc>2 && tdc<12))
1034 //----------------------------------------------------------------
1035 /*void AliTOFQADataMakerRec::ResetAllTRMcounters()
1037 //resets all TRM counters to 0
1038 for (Int_t trm=0;trm<720;trm++){
1039 fTRMNoisyArray[trm]=-1;
1040 fTRMHwOkArray[trm]=-1;
1041 fTRMEnabledArray[trm]=-1;