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),
95 for (Int_t sm=0;sm<10;sm++){
97 fLineSMid3671[sm]=0x0;
102 //____________________________________________________________________________
103 AliTOFQADataMakerRec::AliTOFQADataMakerRec(const AliTOFQADataMakerRec& qadm) :
105 fCalibData(qadm.fCalibData),
106 fEnableNoiseFiltering(qadm.fEnableNoiseFiltering),
107 fEnableDqmShifterOpt(qadm.fEnableDqmShifterOpt),
108 fProcessedRawEventN(qadm.fProcessedRawEventN),
109 fLineExpTimeMin(qadm.fLineExpTimeMin),
110 fLineExpTimeMax(qadm.fLineExpTimeMax),
111 fLineExpTotMin(qadm.fLineExpTotMin),
112 fLineExpTotMax(qadm.fLineExpTotMax)
117 SetName((const char*)qadm.GetName()) ;
118 SetTitle((const char*)qadm.GetTitle());
120 for (Int_t sm=0;sm<10;sm++){
121 fLineSMid035[sm]=qadm.fLineSMid035[sm];
122 fLineSMid3671[sm]=qadm.fLineSMid3671[sm];
126 //__________________________________________________________________
127 AliTOFQADataMakerRec& AliTOFQADataMakerRec::operator = (const AliTOFQADataMakerRec& qadm )
130 // assignment operator.
132 this->~AliTOFQADataMakerRec();
133 new(this) AliTOFQADataMakerRec(qadm);
137 //----------------------------------------------------------------------------
138 AliTOFChannelOnlineStatusArray* AliTOFQADataMakerRec::GetCalibData() const
141 // Retrive TOF calib objects from OCDB
143 AliCDBManager *man = AliCDBManager::Instance();
147 cdbe = man->Get("TOF/Calib/Status",fRun);
149 AliWarning("Load of calibration data from default storage failed!");
150 AliWarning("Calibration data will be loaded from local storage ($ALICE_ROOT)");
151 man->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
152 cdbe = man->Get("TOF/Calib/Status",fRun);
154 // Retrieval of data in directory TOF/Calib/Data:
156 AliTOFChannelOnlineStatusArray * array = 0;
157 if (cdbe) array = (AliTOFChannelOnlineStatusArray *)cdbe->GetObject();
158 if (!array) AliFatal("No calibration data from calibration database !");
166 //____________________________________________________________________________
167 void AliTOFQADataMakerRec::InitRaws()
170 // create Raws histograms in Raws subdir
172 const Bool_t expert = kTRUE ;
173 const Bool_t saveCorr = kTRUE ;
174 const Bool_t image = kTRUE ;
176 TH1I * h0 = new TH1I("hTOFRaws", "TOF raw hit multiplicity; TOF raw hits number;Counts ",200, 0, 200) ;
178 TH1I * h1 = new TH1I("hTOFRawsIA", "TOF raw hit multiplicity - I/A side; TOF raw hits number ;Counts ",200, 0, 200) ;
179 TH1I * h2 = new TH1I("hTOFRawsOA", "TOF raw hit multiplicity - O/A side; TOF raw hits number ;Counts ",200, 0, 200) ;
180 TH1I * h3 = new TH1I("hTOFRawsIC", "TOF raw hit multiplicity - I/C side; TOF raw hits number ;Counts ",200, 0, 200) ;
181 TH1I * h4 = new TH1I("hTOFRawsOC", "TOF raw hit multiplicity - O/C side; TOF raw hits number ;Counts ",200, 0, 200) ;
183 TH1F * h5 = new TH1F("hTOFRawsTime", "TOF Raws - Hit time (ns);Measured Hit time [ns];Counts", 25000,0. ,610.) ;
185 TH1F * h6 = new TH1F("hTOFRawsTimeIA", "TOF Raws - Hit time (ns) - I/A side;Measured Hit time [ns];Counts", 25000,0. ,610.) ;
186 TH1F * h7 = new TH1F("hTOFRawsTimeOA", "TOF Raws - Hit time (ns) - O/A side;Measured Hit time [ns];Counts", 25000,0. ,610.) ;
187 TH1F * h8 = new TH1F("hTOFRawsTimeIC", "TOF Raws - Hit time (ns) - I/C side;Measured Hit time [ns];Counts", 25000,0. ,610.) ;
188 TH1F * h9 = new TH1F("hTOFRawsTimeOC", "TOF Raws - Hit time (ns) - O/C side;Measured Hit time [ns];Counts", 25000,0. ,610.) ;
190 TH1F * h10 = new TH1F("hTOFRawsToT", "TOF Raws - Hit ToT (ns);Measured Hit ToT (ns);Counts", 1000, 0., 48.8) ;
192 TH1F * h11 = new TH1F("hTOFRawsToTIA", "TOF Raws - Hit ToT (ns) - I/A side;Measured Hit ToT (ns);Counts", 1000, 0., 48.8) ;
193 TH1F * h12 = new TH1F("hTOFRawsToTOA", "TOF Raws - Hit ToT (ns) - O/A side;Measured Hit ToT (ns);Counts", 1000, 0., 48.8) ;
194 TH1F * h13 = new TH1F("hTOFRawsToTIC", "TOF Raws - Hit ToT (ns) - I/C side;Measured Hit ToT (ns);Counts", 1000, 0., 48.8) ;
195 TH1F * h14 = new TH1F("hTOFRawsToTOC", "TOF Raws - Hit ToT (ns) - O/C side;Measured Hit ToT (ns);Counts", 1000, 0., 48.8) ;
197 TH1F * h15 = new TH1F("hTOFRawsLTMHits", "LTMs OR signals; Crate; Counts", 72, 0., 72.);
198 TH1F * h16 = new TH1F("hTOFRawsTRMHits035", "TRM hits - crates 0 to 35 ;TRM index = SMid(crate*10)+TRM(0-9);Hits", 361, 0., 361.) ;
199 TH1F * h17 = new TH1F("hTOFRawsTRMHits3671","TRM hits - crates 36 to 71 ;TRM index = SMid(crate*10)+TRM(0-9);Hits", 361, 360., 721.) ;
201 TH1F * h18 = new TH1F("hTOFRawChannelHits","TOF channel hits count; Channel ID; Hits", 158000, 0., 158000);
203 TH1F * h19 = new TH1F("hTOFOrphansTime", "TOF Raws - Orphans time (ns);Measured Hit time [ns];Counts", 25000, 0. ,610.) ;
204 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) ;
205 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) ;
206 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) ;
207 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) ;
208 TH2F * h24 = new TH2F("hTOFTimeVsStrip","TOF Raws - Hit time vs. strip (theta); Strip index;Raws TOF time (ns) ", 91,0.,91, 250, 0., 610.) ;
210 // 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) ;
238 Add2RawsList(h0, 0, !expert, image, !saveCorr) ;
239 Add2RawsList(h1, 1, expert, image, !saveCorr) ;
240 Add2RawsList(h2, 2, expert, image, !saveCorr) ;
241 Add2RawsList(h3, 3, expert, image, !saveCorr) ;
242 Add2RawsList(h4, 4, expert, image, !saveCorr) ;
243 Add2RawsList(h5, 5, !expert, image, !saveCorr) ;
244 Add2RawsList(h6, 6, expert, image, !saveCorr) ;
245 Add2RawsList(h7, 7, expert, image, !saveCorr) ;
246 Add2RawsList(h8, 8, expert, image, !saveCorr) ;
247 Add2RawsList(h9, 9, expert, image, !saveCorr) ;
248 Add2RawsList(h10, 10, !expert, image, !saveCorr) ;
249 Add2RawsList(h11, 11, expert, !image, !saveCorr) ;
250 Add2RawsList(h12, 12, expert, !image, !saveCorr) ;
251 Add2RawsList(h13, 13, expert, !image, !saveCorr) ;
252 Add2RawsList(h14, 14, expert, !image, !saveCorr) ;
253 Add2RawsList(h15, 15, !expert, image, !saveCorr) ;
254 Add2RawsList(h16, 16, !expert, image, !saveCorr) ;
255 Add2RawsList(h17, 17, !expert, image, !saveCorr) ;
256 Add2RawsList(h18, 18, expert, !image, !saveCorr) ;
257 Add2RawsList(h19, 19, expert, !image, !saveCorr) ;
258 Add2RawsList(h20, 20, expert, !image, !saveCorr) ;
259 Add2RawsList(h21, 21, expert, !image, !saveCorr) ;
260 Add2RawsList(h22, 22, expert, !image, !saveCorr) ;
261 Add2RawsList(h23, 23, expert, !image, !saveCorr) ;
262 Add2RawsList(h24, 24, expert, !image, !saveCorr) ;
264 //add lines for DQM shifter
265 fLineExpTimeMin = new TLine(200., 0., 200., 0.);
266 fLineExpTimeMin->SetLineColor(kGreen);
267 fLineExpTimeMin->SetLineWidth(2);
269 fLineExpTimeMax = new TLine(250., 0., 250., 0.);
270 fLineExpTimeMax->SetLineColor(kGreen);
271 fLineExpTimeMax->SetLineWidth(2);
273 fLineExpTotMin = new TLine( 5., 0., 5., 0.);
274 fLineExpTotMin->SetLineColor(kGreen);
275 fLineExpTotMin->SetLineWidth(2);
277 fLineExpTotMax = new TLine(20., 0., 20., 0.);
278 fLineExpTotMax->SetLineColor(kGreen);
279 fLineExpTotMax->SetLineWidth(2);
281 h5->GetListOfFunctions()->Add(fLineExpTimeMin);
282 h5->GetListOfFunctions()->Add(fLineExpTimeMax);
283 h10->GetListOfFunctions()->Add(fLineExpTotMin);
284 h10->GetListOfFunctions()->Add(fLineExpTotMax);
286 for (Int_t sm=0;sm<10;sm++){
287 fLineSMid035[sm] = new TLine( 40*sm, 0, 40*sm, 0.);
288 fLineSMid035[sm]->SetLineColor(kMagenta);
289 fLineSMid035[sm]->SetLineWidth(2);
290 GetRawsData(16)->GetListOfFunctions()->Add(fLineSMid035[sm]);
291 fLineSMid3671[sm] = new TLine( 40*sm+360, 0, 40*sm+360, 0.);
292 fLineSMid3671[sm]->SetLineColor(kMagenta);
293 fLineSMid3671[sm]->SetLineWidth(2);
294 GetRawsData(17)->GetListOfFunctions()->Add(fLineSMid3671[sm]);
299 //____________________________________________________________________________
300 void AliTOFQADataMakerRec::InitRecPoints()
303 // create RecPoints histograms in RecPoints subdir
306 const Bool_t expert = kTRUE ;
307 const Bool_t image = kTRUE ;
309 TH1I * h0 = new TH1I("hTOFRecPoints", "TOF RecPoints multiplicity ; TOF RecPoints number;Counts",200, 0, 200) ;
311 TH1F * h1 = new TH1F("hTOFRecPointsTimeIA", "RecPoints Time Spectrum in TOF (ns)- I/A side; Calibrated TOF time [ns];Counts", 25000, 0., 610.) ;
312 TH1F * h2 = new TH1F("hTOFRecPointsTimeOA", "RecPoints Time Spectrum in TOF (ns)- O/A side; Calibrated TOF time [ns];Counts", 25000, 0., 610.) ;
313 TH1F * h3 = new TH1F("hTOFRecPointsTimeIC", "RecPoints Time Spectrum in TOF (ns)- I/C side; Calibrated TOF time [ns];Counts", 25000, 0., 610.) ;
314 TH1F * h4 = new TH1F("hTOFRecPointsTimeOC", "RecPoints Time Spectrum in TOF (ns)- O/C side; Calibrated TOF time [ns];Counts", 25000, 0., 610.) ;
316 TH1F * h5 = new TH1F("hTOFRecPointsRawTimeIA", "RecPoints raw Time Spectrum in TOF (ns)-I/A side; Measured TOF time [ns];Counts", 25000, 0., 610.) ;
317 TH1F * h6 = new TH1F("hTOFRecPointsRawTimeOA", "RecPoints raw Time Spectrum in TOF (ns)-O/A side; Measured TOF time [ns];Counts", 25000, 0., 610.) ;
318 TH1F * h7 = new TH1F("hTOFRecPointsRawTimeIC", "RecPoints raw Time Spectrum in TOF (ns)-I/C side; Measured TOF time [ns];Counts", 25000, 0., 610.) ;
319 TH1F * h8 = new TH1F("hTOFRecPointsRawTimeOC", "RecPoints raw Time Spectrum in TOF (ns)-O/C side; Measured TOF time [ns];Counts", 25000, 0., 610.) ;
321 TH1F * h9 = new TH1F("hTOFRecPointsToTIA", "RecPoints ToT Spectrum in TOF (ns)-I/A side; Measured TOT [ns];Counts", 1000, 0., 48.8 ) ;
322 TH1F * h10 = new TH1F("hTOFRecPointsToTOA", "RecPoints ToT Spectrum in TOF (ns)-O/A side; Measured TOT [ns];Counts", 1000, 0., 48.8 ) ;
323 TH1F * h11 = new TH1F("hTOFRecPointsToTIC", "RecPoints ToT Spectrum in TOF (ns)-I/C side; Measured TOT [ns];Counts", 1000, 0., 48.8 ) ;
324 TH1F * h12 = new TH1F("hTOFRecPointsToTOC", "RecPoints ToT Spectrum in TOF (ns)-O/C side; Measured TOT [ns];Counts", 1000, 0., 48.8 ) ;
325 TH1F * h13 = new TH1F("hTOFRawChannelHits","TOF channel hits count; Channel ID; Hits", 158000, 0., 158000);
327 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) ;
328 TH2F * h15 = new TH2F("hTOFRecTimeVsStrip","RecPoints TOF time vs. strip (theta); Strip index; RecPoints TOF time (ns) ",92,-1.,91, 250, 0., 610.) ;
330 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) ;
331 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) ;
332 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) ;
333 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) ;
356 Add2RecPointsList(h0, 0, !expert, image) ;
357 Add2RecPointsList(h1, 1, !expert, image) ;
358 Add2RecPointsList(h2, 2, !expert, image) ;
359 Add2RecPointsList(h3, 3, !expert, image) ;
360 Add2RecPointsList(h4, 4, !expert, image) ;
361 Add2RecPointsList(h5, 5, expert, !image) ;
362 Add2RecPointsList(h6, 6, expert, !image) ;
363 Add2RecPointsList(h7, 7, expert, !image) ;
364 Add2RecPointsList(h8, 8, expert, !image) ;
365 Add2RecPointsList(h9, 9, !expert, image) ;
366 Add2RecPointsList(h10, 10, !expert, image) ;
367 Add2RecPointsList(h11, 11, !expert, image) ;
368 Add2RecPointsList(h12, 12, !expert, image) ;
369 Add2RecPointsList(h13, 13, expert, !image) ;
370 Add2RecPointsList(h14, 14, expert, !image) ;
371 Add2RecPointsList(h15, 15, expert, !image) ;
372 Add2RecPointsList(h16, 16, expert, !image) ;
373 Add2RecPointsList(h17, 17, expert, !image) ;
374 Add2RecPointsList(h18, 18, expert, !image) ;
375 Add2RecPointsList(h19, 19, expert, !image) ;
378 //____________________________________________________________________________
379 void AliTOFQADataMakerRec::InitESDs()
382 //create ESDs histograms in ESDs subdir
385 const Bool_t expert = kTRUE ;
386 const Bool_t image = kTRUE ;
388 TH1I * h0 = new TH1I("hTOFESDs", "Number of matched TOF tracks per event;Number of TOF matched ESD tracks;Counts", 200, 0, 200) ;
389 TH1F * h1 = new TH1F("hTOFESDsTime", "Matched ESDs tracks: TOF Time spectrum; Calibrated TOF time [ns];Counts", 25000, 0., 610. ) ;
390 TH1F * h2 = new TH1F("hTOFESDsRawTime", "Matched ESDs tracks: TOF raw Time spectrum;Measured TOF time [ns];Counts", 25000, 0., 610.) ;
391 TH1F * h3 = new TH1F("hTOFESDsToT", "Matched ESDs tracks: TOF ToT spectrum; ESDs ToT [ns];Counts",1000, 0., 48.8) ;
392 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.) ;
393 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.) ;
394 TH1F * h6 = new TH1F("hTOFESDsTPCmatched", "TPC-TOF matched tracks momentum distribution (GeV/c); p (GeV/c);Counts", 50, 0.20, 5.00) ;
395 TH1F * h7 = new TH1F("hTOFESDsMatchingProb", "TPC-TOF track-matching probability per event;TPC-TOF track-matching probability (%) ;Counts",101, -1.0, 100) ;
396 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) ;
397 TH1F * h9 = new TH1F("hTOFESDsDiffTime", "ESDs t_{TOF}-t_{exp} spectrum in TOF (ps); t_{TOF}-t_{exp} [ps];Counts", 200, -2440., 2440.) ;
398 TH1F * h10 = new TH1F("hTOFHitsLength", "Matched ESDs tracks: Length Spectrum; Track length [cm];Counts", 800, 0., 800) ;
411 Add2ESDsList(h0, 0, !expert, image) ;
412 Add2ESDsList(h1, 1, !expert, image) ;
413 Add2ESDsList(h2, 2, expert, image) ;
414 Add2ESDsList(h3, 3, !expert, image) ;
415 Add2ESDsList(h4, 4, expert, image) ;
416 Add2ESDsList(h5, 5, expert, image) ;
417 Add2ESDsList(h6, 6, expert, image) ;
418 Add2ESDsList(h7, 7, expert, image) ;
419 Add2ESDsList(h8, 8, expert, image) ;
420 Add2ESDsList(h9, 9, !expert, image) ;
421 Add2ESDsList(h10, 10, !expert, image) ;
425 //____________________________________________________________________________
426 void AliTOFQADataMakerRec::MakeRaws(AliRawReader* rawReader)
429 // makes data from Raws
431 if (rawReader->GetType()==7) {
433 Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;//in ns
434 Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
436 Int_t ntof[5]; /* 0=tot, 1=IA, 2=OA, 3=IC, 4=OC*/
437 for (Int_t j=0;j<5;j++){ ntof[j]=0;}
439 Int_t equipmentID[5]; //(ddl, trm, chain,tdc,channel)
440 Int_t volumeID[5]; //(sector,plate,strip,padX,padZ)
441 Int_t volumeID2[5]; //(sector,plate,strip,padZ,padX) to use AliTOFGeometry::GetIndex()
442 Int_t out[5]; // out=(indexZ,indexPhi)
445 TClonesArray * clonesRawData;
446 AliTOFRawStream tofInput(rawReader);
448 //uncomment if needed to apply DeltaBC correction
449 //tofInput.ApplyBCCorrections(kTRUE);
451 for (Int_t iDDL = 0; iDDL < AliTOFGeometry::NDDL()*AliTOFGeometry::NSectors(); iDDL++){
454 tofInput.LoadRawDataBuffersV2(iDDL);
455 clonesRawData = (TClonesArray*)tofInput.GetRawData();
456 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
457 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
459 if (tofRawDatum->GetTOF()){
462 equipmentID[1]=tofRawDatum->GetTRM();
463 equipmentID[2]=tofRawDatum->GetTRMchain();
464 equipmentID[3]=tofRawDatum->GetTDC();
465 equipmentID[4]=tofRawDatum->GetTDCchannel();
467 if (CheckEquipID(equipmentID)){
468 tofInput.EquipmentId2VolumeId(iDDL,
469 tofRawDatum->GetTRM(),
470 tofRawDatum->GetTRMchain(),
471 tofRawDatum->GetTDC(),
472 tofRawDatum->GetTDCchannel(),
474 if (FilterSpare(equipmentID)) continue;
475 if (FilterLTMData(equipmentID)){ //counts LTM hits
476 if (tofRawDatum->GetTOT()) GetRawsData(15)->Fill(equipmentID[0]);
478 if (CheckVolumeID(volumeID)){
480 GetMapIndeces(volumeID,out);
481 volumeID2[0]=volumeID[0];
482 volumeID2[1]=volumeID[1];
483 volumeID2[2]=volumeID[2];
484 volumeID2[3]=volumeID[4];
485 volumeID2[4]=volumeID[3];
486 chIndex=AliTOFGeometry::GetIndex(volumeID2);
488 if (tofRawDatum->GetTOT()){
489 if (!(fCalibData->GetNoiseStatus(chIndex)==AliTOFChannelOnlineStatusArray::kTOFNoiseBad)
490 && (fCalibData->GetHWStatus(chIndex) == AliTOFChannelOnlineStatusArray::kTOFHWOk)) {//noise and enabled filter
491 ntof[0]++; //counter for tof hits
493 //fill global spectra for DQM plots
494 GetRawsData(5)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;//in ns
495 GetRawsData(10)->Fill( tofRawDatum->GetTOT()*tot2ns) ;//in ns
497 //fill side-related spectra for experts plots
498 if (volumeID2[0]>4 && volumeID2[0]<14){ //I side
499 if ((iDDL%4==0)|| (iDDL%4==1)){ //A side
501 GetRawsData(6)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;
502 GetRawsData(11)->Fill( tofRawDatum->GetTOT()*tot2ns) ;
504 if ((iDDL%4==2)|| (iDDL%4==3)){//C side
506 GetRawsData(8)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;
507 GetRawsData(13)->Fill( tofRawDatum->GetTOT()*tot2ns) ;
511 if (volumeID2[0]<5 || volumeID2[0]>13){ //O side
512 if ((iDDL%4==0)|| (iDDL%4==1)){ //A side
514 GetRawsData(7)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;
515 GetRawsData(12)->Fill( tofRawDatum->GetTOT()*tot2ns) ;
517 if ((iDDL%4==2)|| (iDDL%4==3)){//C side
519 GetRawsData(9)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;
520 GetRawsData(14)->Fill( tofRawDatum->GetTOT()*tot2ns) ;
525 GetRawsData(18)->Fill(chIndex);
527 Int_t trm= iDDL*10+(equipmentID[1]-3);
528 if (iDDL>=0 && iDDL<36) {
529 GetRawsData(16)->Fill(trm) ;
530 GetRawsData(20)->Fill(trm,tofRawDatum->GetTOF()*tdc2ns);
531 GetRawsData(22)->Fill(trm,tofRawDatum->GetTOT()*tot2ns);
533 if (iDDL>=36 && iDDL<72) {
534 GetRawsData(17)->Fill(trm) ;//in ns
535 GetRawsData(21)->Fill(trm,tofRawDatum->GetTOF()*tdc2ns);
536 GetRawsData(23)->Fill(trm,tofRawDatum->GetTOT()*tot2ns);
538 GetRawsData(24)->Fill(GetStripIndex(volumeID),tofRawDatum->GetTOF()*tdc2ns) ;
539 //GetRawsData(25)->Fill( out[0],out[1]) ;//raw map
543 if (!(fCalibData->GetNoiseStatus(chIndex) == AliTOFChannelOnlineStatusArray::kTOFNoiseBad)
544 && (fCalibData->GetHWStatus(chIndex) == AliTOFChannelOnlineStatusArray::kTOFHWOk))
545 GetRawsData(19)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;//in ns
547 }//end volumeID check
552 clonesRawData->Clear();
555 for (Int_t j=0;j<5;j++){
556 GetRawsData(j)->Fill(ntof[j]);
558 fProcessedRawEventN++;
561 AliDebug(1,Form("Event of type %d found. Skipping non-physics event for QA.\n", rawReader->GetType()));
563 EnableDqmShifterOpt(kTRUE);
566 //____________________________________________________________________________
567 void AliTOFQADataMakerRec::MakeRecPoints(TTree * clustersTree)
570 // Make data from Clusters
573 Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
574 Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
580 TBranch *branch=clustersTree->GetBranch("TOF");
582 AliError("can't get the branch with the TOF clusters !");
586 static TClonesArray dummy("AliTOFcluster",10000);
588 TClonesArray *clusters=&dummy;
589 branch->SetAddress(&clusters);
592 clustersTree->GetEvent(0);
594 GetRecPointsData(0)->Fill((Int_t)clusters->GetEntriesFast()) ;
596 TIter next(clusters) ;
598 while ( (c = dynamic_cast<AliTOFcluster *>(next())) ) {
600 volumeID[0] = c->GetDetInd(0);
601 volumeID[1] = c->GetDetInd(1);
602 volumeID[2] = c->GetDetInd(2);
603 volumeID[3] = c->GetDetInd(4); //X and Z indeces inverted in RecPoints
604 volumeID[4] = c->GetDetInd(3); //X and Z indeces inverted in RecPoints
606 for (Int_t i=0;i<5;i++){
607 volumeID2[i]=c->GetDetInd(i); //X and Z indeces inverted in RecPoints
610 GetMapIndeces(volumeID,out);
611 Int_t chIndex=AliTOFGeometry::GetIndex(volumeID2);
612 Int_t iDDL=AliTOFRawStream::Geant2DDL(volumeID2);
613 Int_t iTRM=AliTOFRawStream::Geant2TRM(volumeID2);
615 if (c->GetTDCRAW() && c->GetTDC() && c->GetToT()){
616 if (volumeID2[0]>4 && volumeID2[0]<14){ //I side
617 if ((iDDL%4==0)|| (iDDL%4==1)){ //A side
618 GetRecPointsData(1)->Fill( c->GetTDC()*tdc2ns) ;//in ns
619 GetRecPointsData(5)->Fill( c->GetTDCRAW()*tdc2ns) ;//in ns
620 GetRecPointsData(9)->Fill( c->GetToT()*tot2ns) ;//in ns
623 if ((iDDL%4==2)|| (iDDL%4==3)){//C side
624 GetRecPointsData(3)->Fill( c->GetTDC()*tdc2ns) ;//in ns
625 GetRecPointsData(7)->Fill( c->GetTDCRAW()*tdc2ns) ;//in ns
626 GetRecPointsData(11)->Fill( c->GetToT()*tot2ns) ;//in ns
630 if (volumeID2[0]<5 || volumeID2[0]>13){ //O side
631 if ((iDDL%4==0)|| (iDDL%4==1)){ //A side
632 GetRecPointsData(2)->Fill( c->GetTDC()*tdc2ns) ;//in ns
633 GetRecPointsData(6)->Fill( c->GetTDCRAW()*tdc2ns) ;//in ns
634 GetRecPointsData(10)->Fill( c->GetToT()*tot2ns) ;//in ns
636 if ((iDDL%4==2)|| (iDDL%4==3)){//C side
637 GetRecPointsData(4)->Fill( c->GetTDC()*tdc2ns) ;//in ns
638 GetRecPointsData(8)->Fill( c->GetTDCRAW()*tdc2ns) ;//in ns
639 GetRecPointsData(12)->Fill( c->GetToT()*tot2ns) ;//in ns
644 GetRecPointsData(14)->Fill(out[0],out[1]);
645 GetRecPointsData(15)->Fill(GetStripIndex(volumeID), c->GetTDC()*tdc2ns) ;
646 Int_t trm= iDDL*10+(iTRM-3);
647 if (iDDL>=0 && iDDL<36) {
648 GetRecPointsData(16)->Fill(trm,c->GetTDC()*tdc2ns);
649 GetRecPointsData(18)->Fill(trm,c->GetToT()*tot2ns);
651 if (iDDL>=36 && iDDL<72) {
652 GetRecPointsData(17)->Fill(trm,c->GetTDC()*tdc2ns);
653 GetRecPointsData(19)->Fill(trm,c->GetToT()*tot2ns);
655 GetRecPointsData(13)->Fill(chIndex) ;//in ns
659 EnableDqmShifterOpt(kFALSE);
662 //____________________________________________________________________________
663 void AliTOFQADataMakerRec::MakeESDs(AliESDEvent * esd)
666 // make QA data from ESDs
668 const Double_t speedOfLight = TMath::C()*1E2*1E-12; // cm/ps
669 const Double_t pionMass = 0.13957018; //GeV/c^2
671 Int_t ntrk = esd->GetNumberOfTracks() ;
680 AliESDtrack *track=esd->GetTrack(ntrk);
681 Double_t tofTime=track->GetTOFsignal();//in ps
682 Double_t tofTimeRaw=track->GetTOFsignalRaw();//in ps
683 Double_t tofToT=track->GetTOFsignalToT(); //in ps
685 UInt_t status=track->GetStatus();
686 if (track->IsOn(AliESDtrack::kTPCrefit)) {
688 Double_t y=track->Eta();
689 if ((status&AliESDtrack::kTOFout)!=0) {
694 if (TMath::Abs(y)<0.9) {
695 GetESDsData(6)->Fill(track->Pt());
698 GetESDsData(1)->Fill(tofTime*1E-3);
699 GetESDsData(2)->Fill(tofTimeRaw*1E-3);
700 GetESDsData(3)->Fill(tofToT*1E-3);
701 //check how many tracks where ESD PID is ok
702 if ((status&AliESDtrack::kESDpid)!=0) nesdpid++;
703 if (((status&AliESDtrack::kESDpid)&AliESDtrack::kTOFpid)!=0) ntofpid++;
705 Double_t length =track->GetIntegratedLength();
706 Double_t mom2=(track->Pt()*track->Pt())+(track->Pz()*track->Pz());
707 Double_t eTexp = TMath::Sqrt(1+(pionMass*pionMass/mom2))*length/speedOfLight; //in ps
708 GetESDsData(9)->Fill(tofTime-eTexp);
709 GetESDsData(10)->Fill(length);
710 } //end check on matched tracks
712 }//end check on TPCrefit
715 GetESDsData(0)->Fill(ntof) ;
719 Float_t ratio = (Int_t)ntofpid/(Int_t)ntof*100; //identified by TOF over matched
720 GetESDsData(4)->Fill(ratio) ;
724 Float_t ratio = (Float_t)ntofpid/(Float_t)nesdpid *100; //identified by TOF over identified
725 GetESDsData(5)->Fill(ratio) ;
729 Float_t ratio = (Float_t)ntof/(Float_t)ntpc*100.; //matching probability
730 GetESDsData(7)->Fill(ratio) ;
734 Float_t ratio = (Float_t)ntof/(Float_t)ntofout*100; //matched over propagated to TOF outer radius
735 GetESDsData(8)->Fill(ratio) ;
737 EnableDqmShifterOpt(kFALSE);
740 //____________________________________________________________________________
741 void AliTOFQADataMakerRec::StartOfDetectorCycle()
744 //Detector specific actions at start of cycle
746 fCalibData = GetCalibData();
750 //____________________________________________________________________________
751 void AliTOFQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
753 //Detector specific actions at end of cycle
754 // do the QA checking
755 for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
756 if ( !AliQAv1::Instance()->IsEventSpecieSet(specie) )
759 if (fEnableDqmShifterOpt){
760 // Help make the raw qa histogram easier to interpret for the DQM shifter
761 if (!GetRawsData(0) || !GetRawsData(5) || !GetRawsData(10)
762 || !GetRawsData(15) || !GetRawsData(16) || !GetRawsData(17)) {
763 printf("No histogram for DQM found - Possible memory corruption ???. Please check\n") ;
766 AliInfo(Form("Processed %i physics raw events",fProcessedRawEventN));
770 //normalize TRM hits plots to the number of enabled channels from OCDB object
771 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.) ;
772 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.) ;
774 for (Int_t ch = 0; ch < fCalibData->GetSize(); ch++) {
775 if (!(fCalibData->GetNoiseStatus(ch)==AliTOFChannelOnlineStatusArray::kTOFNoiseBad)
776 && (fCalibData->GetHWStatus(ch) == AliTOFChannelOnlineStatusArray::kTOFHWOk)){
779 AliTOFGeometry::GetVolumeIndices(ch,geoId);
780 AliTOFRawStream::Geant2EquipmentId(geoId,detId); //detID=(ddl,trm,tdc, chain,channel)
781 if (detId[0]<36) hTrmChannels035->Fill((detId[1]-3)+detId[0]*10);
782 else hTrmChannels3671->Fill((detId[1]-3)+detId[0]*10);
785 GetRawsData(16)->Divide(hTrmChannels035);
786 GetRawsData(16)->SetTitle("TRMs average hit number per active channel - crates 0-35");
787 GetRawsData(16)->GetYaxis()->SetTitle("hits/active channels");
788 GetRawsData(17)->Divide(hTrmChannels3671);
789 GetRawsData(17)->SetTitle("TRMs average hit number per active channel - crates 36-71");
790 GetRawsData(17)->GetYaxis()->SetTitle("hits/active channels");
793 //set minima and maxima to allow log scale
794 Double_t yTimeMax = GetRawsData(5)->GetMaximum();
795 Double_t yTotMax = GetRawsData(10)->GetMaximum();
796 Double_t yTrmMax = TMath::Max(GetRawsData(16)->GetMaximum(),GetRawsData(17)->GetMaximum());
797 // Double_t yLtmMax = GetRawsData(15)->GetMaximum();
798 // Double_t yHitMax = GetRawsData(0)->GetMaximum();
800 // GetRawsData(0)->SetMinimum(0.1);
801 // GetRawsData(5)->SetMinimum(0.1);
802 // GetRawsData(10)->SetMinimum(0.1);
803 // GetRawsData(15)->SetMinimum(0.1);
804 GetRawsData(16)->SetMinimum(0.0001);
805 GetRawsData(17)->SetMinimum(0.0001);
807 // GetRawsData(0)->SetMaximum(yHitMax);
808 // GetRawsData(5)->SetMaximum(yTimeMax);
809 // GetRawsData(10)->SetMaximum(yTotMax);
810 // GetRawsData(15)->SetMaximum(yLtmMax);
811 GetRawsData(16)->SetMaximum(yTrmMax);
812 GetRawsData(17)->SetMaximum(yTrmMax);
814 fLineExpTimeMin->SetY2(yTimeMax);
815 fLineExpTimeMax->SetY2(yTimeMax);
816 fLineExpTotMin->SetY2(yTotMax);
817 fLineExpTotMax->SetY2(yTotMax);
818 for (Int_t sm=0;sm<10;sm++){
819 fLineSMid035[sm]->SetY2(yTrmMax);
820 fLineSMid3671[sm]->SetY2(yTrmMax);
824 //make up for all histos
825 for(Int_t j=0;j<20;j++){
827 GetRawsData(j)->SetMarkerColor(kRed);
828 GetRawsData(j)->SetMarkerStyle(7);
830 GetRawsData(j)->SetLineColor(kBlue+1);
831 GetRawsData(j)->SetLineWidth(1);
832 GetRawsData(j)->SetMarkerColor(kBlue+1);
836 for (Int_t j=15;j<19;j++){
837 GetRawsData(j)->SetFillColor(kGray+1);
838 GetRawsData(j)->SetOption("bar");
843 AliQAChecker::Instance()->Run(AliQAv1::kTOF, task, list) ;
845 //____________________________________________________________________________
846 void AliTOFQADataMakerRec::GetMapIndeces(const Int_t* const in , Int_t* out)
849 //return appropriate indeces for the theta-phi map
852 Int_t npadX = AliTOFGeometry::NpadX();
853 Int_t npadZ = AliTOFGeometry::NpadZ();
854 Int_t nStripA = AliTOFGeometry::NStripA();
855 Int_t nStripB = AliTOFGeometry::NStripB();
856 Int_t nStripC = AliTOFGeometry::NStripC();
858 Int_t isector = in[0];
859 Int_t iplate = in[1];
860 Int_t istrip = in[2];
864 Int_t stripOffset = 0;
870 stripOffset = nStripC;
873 stripOffset = nStripC+nStripB;
876 stripOffset = nStripC+nStripB+nStripA;
879 stripOffset = nStripC+nStripB+nStripA+nStripB;
882 AliDebug(1,Form("Wrong plate number in TOF (%d) !",iplate));
885 Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
886 Int_t phiindex=npadX*isector+ipadX+1;
892 //---------------------------------------------------------------
893 Int_t AliTOFQADataMakerRec::GetStripIndex(const Int_t * const in)
895 /* return tof strip index between 0 and 91 */
897 Int_t nStripA = AliTOFGeometry::NStripA();
898 Int_t nStripB = AliTOFGeometry::NStripB();
899 Int_t nStripC = AliTOFGeometry::NStripC();
901 // Int_t isector = in[0];
902 Int_t iplate = in[1];
903 Int_t istrip = in[2];
904 //Int_t ipadX = in[3];
905 //Int_t ipadZ = in[4];
907 Int_t stripOffset = 0;
913 stripOffset = nStripC;
916 stripOffset = nStripC+nStripB;
919 stripOffset = nStripC+nStripB+nStripA;
922 stripOffset = nStripC+nStripB+nStripA+nStripB;
925 AliDebug(1,Form("Wrong plate number in TOF (%d) !",iplate));
930 if (stripOffset<0 || stripOffset>92) return -1;
932 return (stripOffset+istrip);
935 //---------------------------------------------------------------
936 Bool_t AliTOFQADataMakerRec::CheckVolumeID(const Int_t * const volumeID)
939 //Checks volume ID validity
942 for (Int_t j=0;j<5;j++){
944 AliDebug(1,Form("Invalid detector volume index for volumeID[%i]",j));
952 //---------------------------------------------------------------
953 Bool_t AliTOFQADataMakerRec::CheckEquipID(const Int_t * const equipmentID)
956 //Checks equipment ID validity
958 for (Int_t j=0;j<5;j++){
959 if (equipmentID[j]<0) {
960 AliDebug(1,Form("Invalid equipment volume index for equipmentID[%i]",j));
966 //---------------------------------------------------------------
967 Bool_t AliTOFQADataMakerRec::FilterLTMData(const Int_t * const equipmentID) const
969 /*It returns kTRUE if data come from LTM.
970 It thus filters trigger-related signals */
973 //if (!CheckEquipID(equipmentID)) return kFALSE;
974 ddl = equipmentID[0];
975 trm = equipmentID[1];
976 tdc = equipmentID[3];
978 if ((ddl%2==1) && (trm==3) && (tdc>11 && tdc<15))
984 //---------------------------------------------------------------
985 Bool_t AliTOFQADataMakerRec::FilterSpare(const Int_t * const equipmentID) const
987 /*It returns kTRUE if data come from spare
989 So far only check on TRM 3 crate left is implemented */
992 //if (!CheckEquipID(equipmentID)) return kFALSE;
993 ddl = equipmentID[0];
994 trm = equipmentID[1];
995 tdc = equipmentID[3];
997 if ((ddl%2==1) && (trm==3) && (tdc>2 && tdc<12))