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 13/09/2010
26 - Set TLines as private members
27 - Set image flag for expert histos
29 Modified by fbellini on 14/06/2010
31 - use LoadRawDataBuffersV2()
33 Modified by fbellini on 10/05/2010
34 - Fixed EndOfDetectorCycle() memory corruption bug
36 Modified by fbellini on 22/04/2010
37 - Added filter for physics events
39 Modified by fbellini on 16/04/2010
40 - Added EnableDqmShifterOpt()
41 - Modified EndOfDetectorCycle() with options for DQM
44 Modified by fbellini on 30/03/2010
45 - Changed raws time histos range to 610ns
46 - Added FilterLTMData() and FilterSpare() methods
47 - Added check on enabled channels for raw data
48 - Updated RecPoints QA
50 Modified by fbellini on 02/03/2010
51 - Fixed raw data decoding methods (use AliTOFRawStream::LoadRawDataBuffer())
52 - Added filter for noisy channels and read map from OCDB
53 - Added GetCalibData() method
54 - Added CheckVolumeID() and CheckEquipID() methods
58 #include <TClonesArray.h>
63 #include "AliCDBManager.h"
64 #include "AliCDBEntry.h"
65 #include "AliESDEvent.h"
66 #include "AliESDtrack.h"
67 #include "AliQAChecker.h"
68 #include "AliRawReader.h"
69 #include "AliTOFRawStream.h"
70 #include "AliTOFcluster.h"
71 #include "AliTOFQADataMakerRec.h"
72 #include "AliTOFrawData.h"
73 #include "AliTOFGeometry.h"
74 #include "AliTOFChannelOnlineStatusArray.h"
77 ClassImp(AliTOFQADataMakerRec)
79 //____________________________________________________________________________
80 AliTOFQADataMakerRec::AliTOFQADataMakerRec() :
81 AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kTOF), "TOF Quality Assurance Data Maker"),
83 fEnableNoiseFiltering(kFALSE),
84 fEnableDqmShifterOpt(kFALSE),
85 fProcessedRawEventN(0),
90 fTOFRawStream(AliTOFRawStream())
96 for (Int_t sm=0;sm<10;sm++){
98 fLineSMid3671[sm]=0x0;
103 //____________________________________________________________________________
104 AliTOFQADataMakerRec::AliTOFQADataMakerRec(const AliTOFQADataMakerRec& qadm) :
106 fCalibData(qadm.fCalibData),
107 fEnableNoiseFiltering(qadm.fEnableNoiseFiltering),
108 fEnableDqmShifterOpt(qadm.fEnableDqmShifterOpt),
109 fProcessedRawEventN(qadm.fProcessedRawEventN),
110 fLineExpTimeMin(qadm.fLineExpTimeMin),
111 fLineExpTimeMax(qadm.fLineExpTimeMax),
112 fLineExpTotMin(qadm.fLineExpTotMin),
113 fLineExpTotMax(qadm.fLineExpTotMax),
114 fTOFRawStream(qadm.fTOFRawStream)
119 SetName((const char*)qadm.GetName()) ;
120 SetTitle((const char*)qadm.GetTitle());
122 for (Int_t sm=0;sm<10;sm++){
123 fLineSMid035[sm]=qadm.fLineSMid035[sm];
124 fLineSMid3671[sm]=qadm.fLineSMid3671[sm];
128 //__________________________________________________________________
129 AliTOFQADataMakerRec& AliTOFQADataMakerRec::operator = (const AliTOFQADataMakerRec& qadm )
132 // assignment operator.
134 this->~AliTOFQADataMakerRec();
135 new(this) AliTOFQADataMakerRec(qadm);
139 //----------------------------------------------------------------------------
140 AliTOFQADataMakerRec::~AliTOFQADataMakerRec()
143 fTOFRawStream.Clear();
146 //----------------------------------------------------------------------------
147 AliTOFChannelOnlineStatusArray* AliTOFQADataMakerRec::GetCalibData() const
150 // Retrive TOF calib objects from OCDB
152 AliCDBManager *man = AliCDBManager::Instance();
156 cdbe = man->Get("TOF/Calib/Status",fRun);
158 AliWarning("Load of calibration data from default storage failed!");
159 AliWarning("Calibration data will be loaded from local storage ($ALICE_ROOT)");
160 man->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
161 cdbe = man->Get("TOF/Calib/Status",fRun);
163 // Retrieval of data in directory TOF/Calib/Data:
165 AliTOFChannelOnlineStatusArray * array = 0;
166 if (cdbe) array = (AliTOFChannelOnlineStatusArray *)cdbe->GetObject();
167 if (!array) AliFatal("No calibration data from calibration database !");
175 //____________________________________________________________________________
176 void AliTOFQADataMakerRec::InitRaws()
179 // create Raws histograms in Raws subdir
181 const Bool_t expert = kTRUE ;
182 const Bool_t saveCorr = kTRUE ;
183 const Bool_t image = kTRUE ;
185 TH1I * h0 = new TH1I("hTOFRaws", "TOF raw hit multiplicity; TOF raw hits number;Counts ",200, 0, 200) ;
187 TH1I * h1 = new TH1I("hTOFRawsIA", "TOF raw hit multiplicity - I/A side; TOF raw hits number ;Counts ",200, 0, 200) ;
188 TH1I * h2 = new TH1I("hTOFRawsOA", "TOF raw hit multiplicity - O/A side; TOF raw hits number ;Counts ",200, 0, 200) ;
189 TH1I * h3 = new TH1I("hTOFRawsIC", "TOF raw hit multiplicity - I/C side; TOF raw hits number ;Counts ",200, 0, 200) ;
190 TH1I * h4 = new TH1I("hTOFRawsOC", "TOF raw hit multiplicity - O/C side; TOF raw hits number ;Counts ",200, 0, 200) ;
192 TH1F * h5 = new TH1F("hTOFRawsTime", "TOF Raws - Hit time (ns);Measured Hit time [ns];Counts", 25000,0. ,610.) ;
194 TH1F * h6 = new TH1F("hTOFRawsTimeIA", "TOF Raws - Hit time (ns) - I/A side;Measured Hit time [ns];Counts", 25000,0. ,610.) ;
195 TH1F * h7 = new TH1F("hTOFRawsTimeOA", "TOF Raws - Hit time (ns) - O/A side;Measured Hit time [ns];Counts", 25000,0. ,610.) ;
196 TH1F * h8 = new TH1F("hTOFRawsTimeIC", "TOF Raws - Hit time (ns) - I/C side;Measured Hit time [ns];Counts", 25000,0. ,610.) ;
197 TH1F * h9 = new TH1F("hTOFRawsTimeOC", "TOF Raws - Hit time (ns) - O/C side;Measured Hit time [ns];Counts", 25000,0. ,610.) ;
199 TH1F * h10 = new TH1F("hTOFRawsToT", "TOF Raws - Hit ToT (ns);Measured Hit ToT (ns);Counts", 1000, 0., 48.8) ;
201 TH1F * h11 = new TH1F("hTOFRawsToTIA", "TOF Raws - Hit ToT (ns) - I/A side;Measured Hit ToT (ns);Counts", 1000, 0., 48.8) ;
202 TH1F * h12 = new TH1F("hTOFRawsToTOA", "TOF Raws - Hit ToT (ns) - O/A side;Measured Hit ToT (ns);Counts", 1000, 0., 48.8) ;
203 TH1F * h13 = new TH1F("hTOFRawsToTIC", "TOF Raws - Hit ToT (ns) - I/C side;Measured Hit ToT (ns);Counts", 1000, 0., 48.8) ;
204 TH1F * h14 = new TH1F("hTOFRawsToTOC", "TOF Raws - Hit ToT (ns) - O/C side;Measured Hit ToT (ns);Counts", 1000, 0., 48.8) ;
206 TH1F * h15 = new TH1F("hTOFRawsLTMHits", "LTMs OR signals; Crate; Counts", 72, 0., 72.);
207 TH1F * h16 = new TH1F("hTOFRawsTRMHits035", "TRM hits - crates 0 to 35 ;TRM index = SMid(crate*10)+TRM(0-9);Hits", 361, 0., 361.) ;
208 TH1F * h17 = new TH1F("hTOFRawsTRMHits3671","TRM hits - crates 36 to 71 ;TRM index = SMid(crate*10)+TRM(0-9);Hits", 361, 360., 721.) ;
210 TH1F * h18 = new TH1F("hTOFRawChannelHits","TOF channel hits count; Channel ID; Hits", 158000, 0., 158000);
212 TH1F * h19 = new TH1F("hTOFOrphansTime", "TOF Raws - Orphans time (ns);Measured Hit time [ns];Counts", 25000, 0. ,610.) ;
213 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) ;
214 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) ;
215 TH2F * h22 = new TH2F("hTOFRawToTVsTRM035", "TOF raws - Hit ToT vs TRM - crates 0 to 35; TRM index = DDL*10+TRM(0-9);TOF raw ToT [ns] ", 361, 0., 361., 100, 0., 48.8) ;
216 TH2F * h23 = new TH2F("hTOFRawToTVsTRM3671", "TOF raws - Hit ToT vs TRM - crates 36 to 72; TRM index = DDL*10+TRM(0-9);TOF raw ToT [ns] ", 361, 360., 721., 100, 0., 48.8) ;
217 TH2F * h24 = new TH2F("hTOFTimeVsStrip","TOF Raws - Hit time vs. strip (theta); Strip index;Raws TOF time (ns) ", 91,0.,91, 250, 0., 610.) ;
219 // TH2F * h25 = new TH2F("hTOFRawsClusMap","Raws vs TOF eta-phi;eta (2*strip+padz);phi (48*sector+padx)",183, -0.5, 182.5,865,-0.5,864.5) ;
247 Add2RawsList(h0, 0, !expert, image, !saveCorr) ;
248 Add2RawsList(h1, 1, expert, image, !saveCorr) ;
249 Add2RawsList(h2, 2, expert, image, !saveCorr) ;
250 Add2RawsList(h3, 3, expert, image, !saveCorr) ;
251 Add2RawsList(h4, 4, expert, image, !saveCorr) ;
252 Add2RawsList(h5, 5, !expert, image, !saveCorr) ;
253 Add2RawsList(h6, 6, expert, image, !saveCorr) ;
254 Add2RawsList(h7, 7, expert, image, !saveCorr) ;
255 Add2RawsList(h8, 8, expert, image, !saveCorr) ;
256 Add2RawsList(h9, 9, expert, image, !saveCorr) ;
257 Add2RawsList(h10, 10, !expert, image, !saveCorr) ;
258 Add2RawsList(h11, 11, expert, !image, !saveCorr) ;
259 Add2RawsList(h12, 12, expert, !image, !saveCorr) ;
260 Add2RawsList(h13, 13, expert, !image, !saveCorr) ;
261 Add2RawsList(h14, 14, expert, !image, !saveCorr) ;
262 Add2RawsList(h15, 15, !expert, image, !saveCorr) ;
263 Add2RawsList(h16, 16, !expert, image, !saveCorr) ;
264 Add2RawsList(h17, 17, !expert, image, !saveCorr) ;
265 Add2RawsList(h18, 18, expert, !image, !saveCorr) ;
266 Add2RawsList(h19, 19, expert, !image, !saveCorr) ;
267 Add2RawsList(h20, 20, expert, !image, !saveCorr) ;
268 Add2RawsList(h21, 21, expert, !image, !saveCorr) ;
269 Add2RawsList(h22, 22, expert, !image, !saveCorr) ;
270 Add2RawsList(h23, 23, expert, !image, !saveCorr) ;
271 Add2RawsList(h24, 24, expert, !image, !saveCorr) ;
273 //add lines for DQM shifter
274 fLineExpTimeMin = new TLine(200., 0., 200., 0.);
275 fLineExpTimeMin->SetLineColor(kGreen);
276 fLineExpTimeMin->SetLineWidth(2);
278 fLineExpTimeMax = new TLine(250., 0., 250., 0.);
279 fLineExpTimeMax->SetLineColor(kGreen);
280 fLineExpTimeMax->SetLineWidth(2);
282 fLineExpTotMin = new TLine( 5., 0., 5., 0.);
283 fLineExpTotMin->SetLineColor(kGreen);
284 fLineExpTotMin->SetLineWidth(2);
286 fLineExpTotMax = new TLine(20., 0., 20., 0.);
287 fLineExpTotMax->SetLineColor(kGreen);
288 fLineExpTotMax->SetLineWidth(2);
290 h5->GetListOfFunctions()->Add(fLineExpTimeMin);
291 h5->GetListOfFunctions()->Add(fLineExpTimeMax);
292 h10->GetListOfFunctions()->Add(fLineExpTotMin);
293 h10->GetListOfFunctions()->Add(fLineExpTotMax);
295 for (Int_t sm=0;sm<10;sm++){
296 fLineSMid035[sm] = new TLine( 40*sm, 0, 40*sm, 0.);
297 fLineSMid035[sm]->SetLineColor(kMagenta);
298 fLineSMid035[sm]->SetLineWidth(2);
299 GetRawsData(16)->GetListOfFunctions()->Add(fLineSMid035[sm]);
300 fLineSMid3671[sm] = new TLine( 40*sm+360, 0, 40*sm+360, 0.);
301 fLineSMid3671[sm]->SetLineColor(kMagenta);
302 fLineSMid3671[sm]->SetLineWidth(2);
303 GetRawsData(17)->GetListOfFunctions()->Add(fLineSMid3671[sm]);
308 //____________________________________________________________________________
309 void AliTOFQADataMakerRec::InitRecPoints()
312 // create RecPoints histograms in RecPoints subdir
315 const Bool_t expert = kTRUE ;
316 const Bool_t image = kTRUE ;
318 TH1I * h0 = new TH1I("hTOFRecPoints", "TOF RecPoints multiplicity ; TOF RecPoints number;Counts",200, 0, 200) ;
320 TH1F * h1 = new TH1F("hTOFRecPointsTimeIA", "RecPoints Time Spectrum in TOF (ns)- I/A side; Calibrated TOF time [ns];Counts", 25000, 0., 610.) ;
321 TH1F * h2 = new TH1F("hTOFRecPointsTimeOA", "RecPoints Time Spectrum in TOF (ns)- O/A side; Calibrated TOF time [ns];Counts", 25000, 0., 610.) ;
322 TH1F * h3 = new TH1F("hTOFRecPointsTimeIC", "RecPoints Time Spectrum in TOF (ns)- I/C side; Calibrated TOF time [ns];Counts", 25000, 0., 610.) ;
323 TH1F * h4 = new TH1F("hTOFRecPointsTimeOC", "RecPoints Time Spectrum in TOF (ns)- O/C side; Calibrated TOF time [ns];Counts", 25000, 0., 610.) ;
325 TH1F * h5 = new TH1F("hTOFRecPointsRawTimeIA", "RecPoints raw Time Spectrum in TOF (ns)-I/A side; Measured TOF time [ns];Counts", 25000, 0., 610.) ;
326 TH1F * h6 = new TH1F("hTOFRecPointsRawTimeOA", "RecPoints raw Time Spectrum in TOF (ns)-O/A side; Measured TOF time [ns];Counts", 25000, 0., 610.) ;
327 TH1F * h7 = new TH1F("hTOFRecPointsRawTimeIC", "RecPoints raw Time Spectrum in TOF (ns)-I/C side; Measured TOF time [ns];Counts", 25000, 0., 610.) ;
328 TH1F * h8 = new TH1F("hTOFRecPointsRawTimeOC", "RecPoints raw Time Spectrum in TOF (ns)-O/C side; Measured TOF time [ns];Counts", 25000, 0., 610.) ;
330 TH1F * h9 = new TH1F("hTOFRecPointsToTIA", "RecPoints ToT Spectrum in TOF (ns)-I/A side; Measured TOT [ns];Counts", 1000, 0., 48.8 ) ;
331 TH1F * h10 = new TH1F("hTOFRecPointsToTOA", "RecPoints ToT Spectrum in TOF (ns)-O/A side; Measured TOT [ns];Counts", 1000, 0., 48.8 ) ;
332 TH1F * h11 = new TH1F("hTOFRecPointsToTIC", "RecPoints ToT Spectrum in TOF (ns)-I/C side; Measured TOT [ns];Counts", 1000, 0., 48.8 ) ;
333 TH1F * h12 = new TH1F("hTOFRecPointsToTOC", "RecPoints ToT Spectrum in TOF (ns)-O/C side; Measured TOT [ns];Counts", 1000, 0., 48.8 ) ;
334 TH1F * h13 = new TH1F("hTOFRawChannelHits","TOF channel hits count; Channel ID; Hits", 158000, 0., 158000);
336 TH2F * h14 = new TH2F("hTOFRecPointsClusMap","RecPoints vs TOF phi-eta; eta (2*strip+padz);phi (48*sector+padx)",183, -0.5, 182.5,865,-0.5,864.5) ;
337 TH2F * h15 = new TH2F("hTOFRecTimeVsStrip","RecPoints TOF time vs. strip (theta); Strip index; RecPoints TOF time (ns) ",92,-1.,91, 250, 0., 610.) ;
339 TH2F * h16 = 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) ;
340 TH2F * h17 = 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) ;
341 TH2F * h18 = new TH2F("hTOFRecPointsToTVsTRM035","TOF RecPoints ToT vs TRM - crates 0 to 35; TRM index = DDL*10+TRM(0-9);TOF ToT [ns] ", 361, 0., 361., 100, 0., 48.8) ;
342 TH2F * h19 = new TH2F("hTOFRecPointsToTVsTRM3671","TOF RecPoints ToT vs TRM - crates 36 to 72; TRM index = DDL*10+TRM(0-9);TOF ToT [ns] ", 361, 360., 721., 100, 0., 48.8) ;
365 Add2RecPointsList(h0, 0, !expert, image) ;
366 Add2RecPointsList(h1, 1, !expert, image) ;
367 Add2RecPointsList(h2, 2, !expert, image) ;
368 Add2RecPointsList(h3, 3, !expert, image) ;
369 Add2RecPointsList(h4, 4, !expert, image) ;
370 Add2RecPointsList(h5, 5, expert, !image) ;
371 Add2RecPointsList(h6, 6, expert, !image) ;
372 Add2RecPointsList(h7, 7, expert, !image) ;
373 Add2RecPointsList(h8, 8, expert, !image) ;
374 Add2RecPointsList(h9, 9, !expert, image) ;
375 Add2RecPointsList(h10, 10, !expert, image) ;
376 Add2RecPointsList(h11, 11, !expert, image) ;
377 Add2RecPointsList(h12, 12, !expert, image) ;
378 Add2RecPointsList(h13, 13, expert, !image) ;
379 Add2RecPointsList(h14, 14, expert, !image) ;
380 Add2RecPointsList(h15, 15, expert, !image) ;
381 Add2RecPointsList(h16, 16, expert, !image) ;
382 Add2RecPointsList(h17, 17, expert, !image) ;
383 Add2RecPointsList(h18, 18, expert, !image) ;
384 Add2RecPointsList(h19, 19, expert, !image) ;
387 //____________________________________________________________________________
388 void AliTOFQADataMakerRec::InitESDs()
391 //create ESDs histograms in ESDs subdir
394 const Bool_t expert = kTRUE ;
395 const Bool_t image = kTRUE ;
397 TH1I * h0 = new TH1I("hTOFESDs", "Number of matched TOF tracks per event;Number of TOF matched ESD tracks;Counts", 200, 0, 200) ;
398 TH1F * h1 = new TH1F("hTOFESDsTime", "Matched ESDs tracks: TOF Time spectrum; Calibrated TOF time [ns];Counts", 25000, 0., 610. ) ;
399 TH1F * h2 = new TH1F("hTOFESDsRawTime", "Matched ESDs tracks: TOF raw Time spectrum;Measured TOF time [ns];Counts", 25000, 0., 610.) ;
400 TH1F * h3 = new TH1F("hTOFESDsToT", "Matched ESDs tracks: TOF ToT spectrum; ESDs ToT [ns];Counts",1000, 0., 48.8) ;
401 TH1F * h4 = new TH1F("hTOFESDsPIDoverMatched", "Fraction of TOF identified tracks over matched TOF tracks per event; identified by TOF / matched tracks [%];Counts", 101, -1., 100.) ;
402 TH1F * h5 = new TH1F("hTOFESDsPID", "Fraction of TOF identified tracks over ESD identified tracks per event; identified by TOF / ESD identified [%];Counts", 101, -1., 100.) ;
403 TH1F * h6 = new TH1F("hTOFESDsTPCmatched", "TPC-TOF matched tracks momentum distribution (GeV/c); p (GeV/c);Counts", 50, 0.20, 5.00) ;
404 TH1F * h7 = new TH1F("hTOFESDsMatchingProb", "TPC-TOF track-matching probability per event;TPC-TOF track-matching probability (%) ;Counts",101, -1.0, 100) ;
405 TH1F * h8 = new TH1F("hTOFESDsHitOverTracked", "Fraction of TOF matching tracks over propagated tracks per event; TOF matched / propagated to TOF [%];Counts",101, -1.0, 100.0) ;
406 TH1F * h9 = new TH1F("hTOFESDsDiffTime", "ESDs t_{TOF}-t_{exp} spectrum in TOF (ps); t_{TOF}-t_{exp} [ps];Counts", 200, -2440., 2440.) ;
407 TH1F * h10 = new TH1F("hTOFHitsLength", "Matched ESDs tracks: Length Spectrum; Track length [cm];Counts", 800, 0., 800) ;
420 Add2ESDsList(h0, 0, !expert, image) ;
421 Add2ESDsList(h1, 1, !expert, image) ;
422 Add2ESDsList(h2, 2, expert, image) ;
423 Add2ESDsList(h3, 3, !expert, image) ;
424 Add2ESDsList(h4, 4, expert, image) ;
425 Add2ESDsList(h5, 5, expert, image) ;
426 Add2ESDsList(h6, 6, expert, image) ;
427 Add2ESDsList(h7, 7, expert, image) ;
428 Add2ESDsList(h8, 8, expert, image) ;
429 Add2ESDsList(h9, 9, !expert, image) ;
430 Add2ESDsList(h10, 10, !expert, image) ;
434 //____________________________________________________________________________
435 void AliTOFQADataMakerRec::MakeRaws(AliRawReader* rawReader)
438 // makes data from Raws
440 if (rawReader->GetType()==7) {
442 Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;//in ns
443 Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
445 Int_t ntof[5]; /* 0=tot, 1=IA, 2=OA, 3=IC, 4=OC*/
446 for (Int_t j=0;j<5;j++){ ntof[j]=0;}
448 Int_t equipmentID[5]; //(ddl, trm, chain,tdc,channel)
449 Int_t volumeID[5]; //(sector,plate,strip,padX,padZ)
450 Int_t volumeID2[5]; //(sector,plate,strip,padZ,padX) to use AliTOFGeometry::GetIndex()
451 Int_t out[5]; // out=(indexZ,indexPhi)
454 TClonesArray * clonesRawData;
455 //AliTOFRawStream tofInput(rawReader);
456 fTOFRawStream.SetRawReader(rawReader);
458 //uncomment if needed to apply DeltaBC correction
459 //fTOFRawStream.ApplyBCCorrections(kTRUE);
461 for (Int_t iDDL = 0; iDDL < AliTOFGeometry::NDDL()*AliTOFGeometry::NSectors(); iDDL++){
464 fTOFRawStream.LoadRawDataBuffersV2(iDDL);
465 clonesRawData = (TClonesArray*)fTOFRawStream.GetRawData();
466 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
467 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
469 if (tofRawDatum->GetTOF()){
472 equipmentID[1]=tofRawDatum->GetTRM();
473 equipmentID[2]=tofRawDatum->GetTRMchain();
474 equipmentID[3]=tofRawDatum->GetTDC();
475 equipmentID[4]=tofRawDatum->GetTDCchannel();
477 if (CheckEquipID(equipmentID)){
478 fTOFRawStream.EquipmentId2VolumeId(iDDL,
479 tofRawDatum->GetTRM(),
480 tofRawDatum->GetTRMchain(),
481 tofRawDatum->GetTDC(),
482 tofRawDatum->GetTDCchannel(),
484 if (FilterSpare(equipmentID)) continue;
485 if (FilterLTMData(equipmentID)){ //counts LTM hits
486 if (tofRawDatum->GetTOT()) GetRawsData(15)->Fill(equipmentID[0]);
488 if (CheckVolumeID(volumeID)){
490 GetMapIndeces(volumeID,out);
491 volumeID2[0]=volumeID[0];
492 volumeID2[1]=volumeID[1];
493 volumeID2[2]=volumeID[2];
494 volumeID2[3]=volumeID[4];
495 volumeID2[4]=volumeID[3];
496 chIndex=AliTOFGeometry::GetIndex(volumeID2);
498 if (tofRawDatum->GetTOT()){
499 if (!(fCalibData->GetNoiseStatus(chIndex)==AliTOFChannelOnlineStatusArray::kTOFNoiseBad)
500 && (fCalibData->GetHWStatus(chIndex) == AliTOFChannelOnlineStatusArray::kTOFHWOk)) {//noise and enabled filter
501 ntof[0]++; //counter for tof hits
503 //fill global spectra for DQM plots
504 GetRawsData(5)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;//in ns
505 GetRawsData(10)->Fill( tofRawDatum->GetTOT()*tot2ns) ;//in ns
507 //fill side-related spectra for experts plots
508 if (volumeID2[0]>4 && volumeID2[0]<14){ //I side
509 if ((iDDL%4==0)|| (iDDL%4==1)){ //A side
511 GetRawsData(6)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;
512 GetRawsData(11)->Fill( tofRawDatum->GetTOT()*tot2ns) ;
514 if ((iDDL%4==2)|| (iDDL%4==3)){//C side
516 GetRawsData(8)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;
517 GetRawsData(13)->Fill( tofRawDatum->GetTOT()*tot2ns) ;
521 if (volumeID2[0]<5 || volumeID2[0]>13){ //O side
522 if ((iDDL%4==0)|| (iDDL%4==1)){ //A side
524 GetRawsData(7)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;
525 GetRawsData(12)->Fill( tofRawDatum->GetTOT()*tot2ns) ;
527 if ((iDDL%4==2)|| (iDDL%4==3)){//C side
529 GetRawsData(9)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;
530 GetRawsData(14)->Fill( tofRawDatum->GetTOT()*tot2ns) ;
535 GetRawsData(18)->Fill(chIndex);
537 Int_t trm= iDDL*10+(equipmentID[1]-3);
538 if (iDDL>=0 && iDDL<36) {
539 GetRawsData(16)->Fill(trm) ;
540 GetRawsData(20)->Fill(trm,tofRawDatum->GetTOF()*tdc2ns);
541 GetRawsData(22)->Fill(trm,tofRawDatum->GetTOT()*tot2ns);
543 if (iDDL>=36 && iDDL<72) {
544 GetRawsData(17)->Fill(trm) ;//in ns
545 GetRawsData(21)->Fill(trm,tofRawDatum->GetTOF()*tdc2ns);
546 GetRawsData(23)->Fill(trm,tofRawDatum->GetTOT()*tot2ns);
548 GetRawsData(24)->Fill(GetStripIndex(volumeID),tofRawDatum->GetTOF()*tdc2ns) ;
549 //GetRawsData(25)->Fill( out[0],out[1]) ;//raw map
553 if (!(fCalibData->GetNoiseStatus(chIndex) == AliTOFChannelOnlineStatusArray::kTOFNoiseBad)
554 && (fCalibData->GetHWStatus(chIndex) == AliTOFChannelOnlineStatusArray::kTOFHWOk))
555 GetRawsData(19)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;//in ns
557 }//end volumeID check
562 clonesRawData->Clear();
565 for (Int_t j=0;j<5;j++){
566 GetRawsData(j)->Fill(ntof[j]);
568 fProcessedRawEventN++;
570 fTOFRawStream.Clear();
573 AliDebug(1,Form("Event of type %d found. Skipping non-physics event for QA.\n", rawReader->GetType()));
575 EnableDqmShifterOpt(kTRUE);
578 //____________________________________________________________________________
579 void AliTOFQADataMakerRec::MakeRecPoints(TTree * clustersTree)
582 // Make data from Clusters
585 Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
586 Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
592 TBranch *branch=clustersTree->GetBranch("TOF");
594 AliError("can't get the branch with the TOF clusters !");
598 static TClonesArray dummy("AliTOFcluster",10000);
600 TClonesArray *clusters=&dummy;
601 branch->SetAddress(&clusters);
604 clustersTree->GetEvent(0);
606 GetRecPointsData(0)->Fill((Int_t)clusters->GetEntriesFast()) ;
608 TIter next(clusters) ;
610 while ( (c = dynamic_cast<AliTOFcluster *>(next())) ) {
612 volumeID[0] = c->GetDetInd(0);
613 volumeID[1] = c->GetDetInd(1);
614 volumeID[2] = c->GetDetInd(2);
615 volumeID[3] = c->GetDetInd(4); //X and Z indeces inverted in RecPoints
616 volumeID[4] = c->GetDetInd(3); //X and Z indeces inverted in RecPoints
618 for (Int_t i=0;i<5;i++){
619 volumeID2[i]=c->GetDetInd(i); //X and Z indeces inverted in RecPoints
622 GetMapIndeces(volumeID,out);
623 Int_t chIndex=AliTOFGeometry::GetIndex(volumeID2);
624 Int_t iDDL=AliTOFRawStream::Geant2DDL(volumeID2);
625 Int_t iTRM=AliTOFRawStream::Geant2TRM(volumeID2);
627 if (c->GetTDCRAW() && c->GetTDC() && c->GetToT()){
628 if (volumeID2[0]>4 && volumeID2[0]<14){ //I side
629 if ((iDDL%4==0)|| (iDDL%4==1)){ //A side
630 GetRecPointsData(1)->Fill( c->GetTDC()*tdc2ns) ;//in ns
631 GetRecPointsData(5)->Fill( c->GetTDCRAW()*tdc2ns) ;//in ns
632 GetRecPointsData(9)->Fill( c->GetToT()*tot2ns) ;//in ns
635 if ((iDDL%4==2)|| (iDDL%4==3)){//C side
636 GetRecPointsData(3)->Fill( c->GetTDC()*tdc2ns) ;//in ns
637 GetRecPointsData(7)->Fill( c->GetTDCRAW()*tdc2ns) ;//in ns
638 GetRecPointsData(11)->Fill( c->GetToT()*tot2ns) ;//in ns
642 if (volumeID2[0]<5 || volumeID2[0]>13){ //O side
643 if ((iDDL%4==0)|| (iDDL%4==1)){ //A side
644 GetRecPointsData(2)->Fill( c->GetTDC()*tdc2ns) ;//in ns
645 GetRecPointsData(6)->Fill( c->GetTDCRAW()*tdc2ns) ;//in ns
646 GetRecPointsData(10)->Fill( c->GetToT()*tot2ns) ;//in ns
648 if ((iDDL%4==2)|| (iDDL%4==3)){//C side
649 GetRecPointsData(4)->Fill( c->GetTDC()*tdc2ns) ;//in ns
650 GetRecPointsData(8)->Fill( c->GetTDCRAW()*tdc2ns) ;//in ns
651 GetRecPointsData(12)->Fill( c->GetToT()*tot2ns) ;//in ns
656 GetRecPointsData(14)->Fill(out[0],out[1]);
657 GetRecPointsData(15)->Fill(GetStripIndex(volumeID), c->GetTDC()*tdc2ns) ;
658 Int_t trm= iDDL*10+(iTRM-3);
659 if (iDDL>=0 && iDDL<36) {
660 GetRecPointsData(16)->Fill(trm,c->GetTDC()*tdc2ns);
661 GetRecPointsData(18)->Fill(trm,c->GetToT()*tot2ns);
663 if (iDDL>=36 && iDDL<72) {
664 GetRecPointsData(17)->Fill(trm,c->GetTDC()*tdc2ns);
665 GetRecPointsData(19)->Fill(trm,c->GetToT()*tot2ns);
667 GetRecPointsData(13)->Fill(chIndex) ;//in ns
671 EnableDqmShifterOpt(kFALSE);
674 //____________________________________________________________________________
675 void AliTOFQADataMakerRec::MakeESDs(AliESDEvent * esd)
678 // make QA data from ESDs
680 const Double_t speedOfLight = TMath::C()*1E2*1E-12; // cm/ps
681 const Double_t pionMass = 0.13957018; //GeV/c^2
683 Int_t ntrk = esd->GetNumberOfTracks() ;
692 AliESDtrack *track=esd->GetTrack(ntrk);
693 Double_t tofTime=track->GetTOFsignal();//in ps
694 Double_t tofTimeRaw=track->GetTOFsignalRaw();//in ps
695 Double_t tofToT=track->GetTOFsignalToT(); //in ps
697 UInt_t status=track->GetStatus();
698 if (track->IsOn(AliESDtrack::kTPCrefit)) {
700 Double_t y=track->Eta();
701 if ((status&AliESDtrack::kTOFout)!=0) {
706 if (TMath::Abs(y)<0.9) {
707 GetESDsData(6)->Fill(track->Pt());
710 GetESDsData(1)->Fill(tofTime*1E-3);
711 GetESDsData(2)->Fill(tofTimeRaw*1E-3);
712 GetESDsData(3)->Fill(tofToT*1E-3);
713 //check how many tracks where ESD PID is ok
714 if ((status&AliESDtrack::kESDpid)!=0) nesdpid++;
715 if (((status&AliESDtrack::kESDpid)&AliESDtrack::kTOFpid)!=0) ntofpid++;
717 Double_t length =track->GetIntegratedLength();
718 Double_t mom2=(track->Pt()*track->Pt())+(track->Pz()*track->Pz());
719 Double_t eTexp = TMath::Sqrt(1+(pionMass*pionMass/mom2))*length/speedOfLight; //in ps
720 GetESDsData(9)->Fill(tofTime-eTexp);
721 GetESDsData(10)->Fill(length);
722 } //end check on matched tracks
724 }//end check on TPCrefit
727 GetESDsData(0)->Fill(ntof) ;
731 Float_t ratio = (Int_t)ntofpid/(Int_t)ntof*100; //identified by TOF over matched
732 GetESDsData(4)->Fill(ratio) ;
736 Float_t ratio = (Float_t)ntofpid/(Float_t)nesdpid *100; //identified by TOF over identified
737 GetESDsData(5)->Fill(ratio) ;
741 Float_t ratio = (Float_t)ntof/(Float_t)ntpc*100.; //matching probability
742 GetESDsData(7)->Fill(ratio) ;
746 Float_t ratio = (Float_t)ntof/(Float_t)ntofout*100; //matched over propagated to TOF outer radius
747 GetESDsData(8)->Fill(ratio) ;
749 EnableDqmShifterOpt(kFALSE);
752 //____________________________________________________________________________
753 void AliTOFQADataMakerRec::StartOfDetectorCycle()
756 //Detector specific actions at start of cycle
758 fCalibData = GetCalibData();
762 //____________________________________________________________________________
763 void AliTOFQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
765 //Detector specific actions at end of cycle
766 // do the QA checking
767 for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
768 if ( !AliQAv1::Instance()->IsEventSpecieSet(specie) )
771 if (fEnableDqmShifterOpt){
772 // Help make the raw qa histogram easier to interpret for the DQM shifter
773 if (!GetRawsData(0) || !GetRawsData(5) || !GetRawsData(10)
774 || !GetRawsData(15) || !GetRawsData(16) || !GetRawsData(17)) {
775 printf("No histogram for DQM found - Possible memory corruption ???. Please check\n") ;
778 AliInfo(Form("Processed %i physics raw events",fProcessedRawEventN));
782 //normalize TRM hits plots to the number of enabled channels from OCDB object
783 TH1F * hTrmChannels035 = new TH1F("hTrmchannels035", "Active channels per TRM - crates 0 to 35;TRM index = SMid(crate*10)+TRM(0-9);Active channels", 361, 0., 361.) ;
784 TH1F * hTrmChannels3671 = new TH1F("hTrmChannels3671","Active channels per TRM - crates 36 to 71 ;TRM index = SMid(crate*10)+TRM(0-9);Active channels", 361, 360., 721.) ;
786 for (Int_t ch = 0; ch < fCalibData->GetSize(); ch++) {
787 if (!(fCalibData->GetNoiseStatus(ch)==AliTOFChannelOnlineStatusArray::kTOFNoiseBad)
788 && (fCalibData->GetHWStatus(ch) == AliTOFChannelOnlineStatusArray::kTOFHWOk)){
791 AliTOFGeometry::GetVolumeIndices(ch,geoId);
792 AliTOFRawStream::Geant2EquipmentId(geoId,detId); //detID=(ddl,trm,tdc, chain,channel)
793 if (detId[0]<36) hTrmChannels035->Fill((detId[1]-3)+detId[0]*10);
794 else hTrmChannels3671->Fill((detId[1]-3)+detId[0]*10);
797 GetRawsData(16)->Divide(hTrmChannels035);
798 GetRawsData(16)->SetTitle("TRMs average hit number per active channel - crates 0-35");
799 GetRawsData(16)->GetYaxis()->SetTitle("hits/active channels");
800 GetRawsData(17)->Divide(hTrmChannels3671);
801 GetRawsData(17)->SetTitle("TRMs average hit number per active channel - crates 36-71");
802 GetRawsData(17)->GetYaxis()->SetTitle("hits/active channels");
805 //set minima and maxima to allow log scale
806 Double_t yTimeMax = GetRawsData(5)->GetMaximum();
807 Double_t yTotMax = GetRawsData(10)->GetMaximum();
808 Double_t yTrmMax = TMath::Max(GetRawsData(16)->GetMaximum(),GetRawsData(17)->GetMaximum());
809 // Double_t yLtmMax = GetRawsData(15)->GetMaximum();
810 // Double_t yHitMax = GetRawsData(0)->GetMaximum();
812 // GetRawsData(0)->SetMinimum(0.1);
813 // GetRawsData(5)->SetMinimum(0.1);
814 // GetRawsData(10)->SetMinimum(0.1);
815 // GetRawsData(15)->SetMinimum(0.1);
816 GetRawsData(16)->SetMinimum(0.0001);
817 GetRawsData(17)->SetMinimum(0.0001);
819 // GetRawsData(0)->SetMaximum(yHitMax);
820 // GetRawsData(5)->SetMaximum(yTimeMax);
821 // GetRawsData(10)->SetMaximum(yTotMax);
822 // GetRawsData(15)->SetMaximum(yLtmMax);
823 GetRawsData(16)->SetMaximum(yTrmMax);
824 GetRawsData(17)->SetMaximum(yTrmMax);
826 fLineExpTimeMin->SetY2(yTimeMax);
827 fLineExpTimeMax->SetY2(yTimeMax);
828 fLineExpTotMin->SetY2(yTotMax);
829 fLineExpTotMax->SetY2(yTotMax);
830 for (Int_t sm=0;sm<10;sm++){
831 fLineSMid035[sm]->SetY2(yTrmMax);
832 fLineSMid3671[sm]->SetY2(yTrmMax);
836 //make up for all histos
837 for(Int_t j=0;j<20;j++){
839 GetRawsData(j)->SetMarkerColor(kRed);
840 GetRawsData(j)->SetMarkerStyle(7);
842 GetRawsData(j)->SetLineColor(kBlue+1);
843 GetRawsData(j)->SetLineWidth(1);
844 GetRawsData(j)->SetMarkerColor(kBlue+1);
848 for (Int_t j=15;j<19;j++){
849 GetRawsData(j)->SetFillColor(kGray+1);
850 GetRawsData(j)->SetOption("bar");
855 AliQAChecker::Instance()->Run(AliQAv1::kTOF, task, list) ;
857 //____________________________________________________________________________
858 void AliTOFQADataMakerRec::GetMapIndeces(const Int_t* const in , Int_t* out)
861 //return appropriate indeces for the theta-phi map
864 Int_t npadX = AliTOFGeometry::NpadX();
865 Int_t npadZ = AliTOFGeometry::NpadZ();
866 Int_t nStripA = AliTOFGeometry::NStripA();
867 Int_t nStripB = AliTOFGeometry::NStripB();
868 Int_t nStripC = AliTOFGeometry::NStripC();
870 Int_t isector = in[0];
871 Int_t iplate = in[1];
872 Int_t istrip = in[2];
876 Int_t stripOffset = 0;
882 stripOffset = nStripC;
885 stripOffset = nStripC+nStripB;
888 stripOffset = nStripC+nStripB+nStripA;
891 stripOffset = nStripC+nStripB+nStripA+nStripB;
894 AliDebug(1,Form("Wrong plate number in TOF (%d) !",iplate));
897 Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
898 Int_t phiindex=npadX*isector+ipadX+1;
904 //---------------------------------------------------------------
905 Int_t AliTOFQADataMakerRec::GetStripIndex(const Int_t * const in)
907 /* return tof strip index between 0 and 91 */
909 Int_t nStripA = AliTOFGeometry::NStripA();
910 Int_t nStripB = AliTOFGeometry::NStripB();
911 Int_t nStripC = AliTOFGeometry::NStripC();
913 // Int_t isector = in[0];
914 Int_t iplate = in[1];
915 Int_t istrip = in[2];
916 //Int_t ipadX = in[3];
917 //Int_t ipadZ = in[4];
919 Int_t stripOffset = 0;
925 stripOffset = nStripC;
928 stripOffset = nStripC+nStripB;
931 stripOffset = nStripC+nStripB+nStripA;
934 stripOffset = nStripC+nStripB+nStripA+nStripB;
937 AliDebug(1,Form("Wrong plate number in TOF (%d) !",iplate));
942 if (stripOffset<0 || stripOffset>92) return -1;
944 return (stripOffset+istrip);
947 //---------------------------------------------------------------
948 Bool_t AliTOFQADataMakerRec::CheckVolumeID(const Int_t * const volumeID)
951 //Checks volume ID validity
954 for (Int_t j=0;j<5;j++){
956 AliDebug(1,Form("Invalid detector volume index for volumeID[%i]",j));
964 //---------------------------------------------------------------
965 Bool_t AliTOFQADataMakerRec::CheckEquipID(const Int_t * const equipmentID)
968 //Checks equipment ID validity
970 for (Int_t j=0;j<5;j++){
971 if (equipmentID[j]<0) {
972 AliDebug(1,Form("Invalid equipment volume index for equipmentID[%i]",j));
978 //---------------------------------------------------------------
979 Bool_t AliTOFQADataMakerRec::FilterLTMData(const Int_t * const equipmentID) const
981 /*It returns kTRUE if data come from LTM.
982 It thus filters trigger-related signals */
985 //if (!CheckEquipID(equipmentID)) return kFALSE;
986 ddl = equipmentID[0];
987 trm = equipmentID[1];
988 tdc = equipmentID[3];
990 if ((ddl%2==1) && (trm==3) && (tdc>11 && tdc<15))
996 //---------------------------------------------------------------
997 Bool_t AliTOFQADataMakerRec::FilterSpare(const Int_t * const equipmentID) const
999 /*It returns kTRUE if data come from spare
1001 So far only check on TRM 3 crate left is implemented */
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>2 && tdc<12))