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 10/05/2010
26 - Fixed EndOfDetectorCycle() memory corruption bug
28 Modified by fbellini on 22/04/2010
29 - Added filter for physics events
31 Modified by fbellini on 16/04/2010
32 - Added EnableDqmShifterOpt()
33 - Modified EndOfDetectorCycle() with options for DQM
36 Modified by fbellini on 30/03/2010
37 - Changed raws time histos range to 610ns
38 - Added FilterLTMData() and FilterSpare() methods
39 - Added check on enabled channels for raw data
40 - Updated RecPoints QA
42 Modified by fbellini on 02/03/2010
43 - Fixed raw data decoding methods (use AliTOFRawStream::LoadRawDataBuffer())
44 - Added filter for noisy channels and read map from OCDB
45 - Added GetCalibData() method
46 - Added CheckVolumeID() and CheckEquipID() methods
50 #include <TClonesArray.h>
55 #include "AliCDBManager.h"
56 #include "AliCDBEntry.h"
57 #include "AliESDEvent.h"
58 #include "AliESDtrack.h"
59 #include "AliQAChecker.h"
60 #include "AliRawReader.h"
61 #include "AliTOFRawStream.h"
62 #include "AliTOFcluster.h"
63 #include "AliTOFQADataMakerRec.h"
64 #include "AliTOFrawData.h"
65 #include "AliTOFGeometry.h"
66 #include "AliTOFChannelOnlineStatusArray.h"
69 ClassImp(AliTOFQADataMakerRec)
71 //____________________________________________________________________________
72 AliTOFQADataMakerRec::AliTOFQADataMakerRec() :
73 AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kTOF), "TOF Quality Assurance Data Maker"),
75 fEnableNoiseFiltering(kFALSE),
76 fEnableDqmShifterOpt(kFALSE),
77 fProcessedRawEventN(0)
85 //____________________________________________________________________________
86 AliTOFQADataMakerRec::AliTOFQADataMakerRec(const AliTOFQADataMakerRec& qadm) :
88 fCalibData(qadm.fCalibData),
89 fEnableNoiseFiltering(qadm.fEnableNoiseFiltering),
90 fEnableDqmShifterOpt(qadm.fEnableDqmShifterOpt),
91 fProcessedRawEventN(qadm.fProcessedRawEventN)
96 SetName((const char*)qadm.GetName()) ;
97 SetTitle((const char*)qadm.GetTitle());
102 //__________________________________________________________________
103 AliTOFQADataMakerRec& AliTOFQADataMakerRec::operator = (const AliTOFQADataMakerRec& qadm )
106 // assignment operator.
108 this->~AliTOFQADataMakerRec();
109 new(this) AliTOFQADataMakerRec(qadm);
113 //----------------------------------------------------------------------------
114 AliTOFChannelOnlineStatusArray* AliTOFQADataMakerRec::GetCalibData() const
117 // Retrive TOF calib objects from OCDB
119 AliCDBManager *man = AliCDBManager::Instance();
123 cdbe = man->Get("TOF/Calib/Status",fRun);
125 AliWarning("Load of calibration data from default storage failed!");
126 AliWarning("Calibration data will be loaded from local storage ($ALICE_ROOT)");
127 man->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
128 cdbe = man->Get("TOF/Calib/Status",fRun);
130 // Retrieval of data in directory TOF/Calib/Data:
132 AliTOFChannelOnlineStatusArray * array = 0;
133 if (cdbe) array = (AliTOFChannelOnlineStatusArray *)cdbe->GetObject();
134 if (!array) AliFatal("No calibration data from calibration database !");
142 //____________________________________________________________________________
143 void AliTOFQADataMakerRec::InitRaws()
146 // create Raws histograms in Raws subdir
148 const Bool_t expert = kTRUE ;
149 const Bool_t saveCorr = kTRUE ;
150 const Bool_t image = kTRUE ;
152 TH1I * h0 = new TH1I("hTOFRaws", "TOF raw hit multiplicity; TOF raw hits number;Counts ",200, 0, 200) ;
154 TH1I * h1 = new TH1I("hTOFRawsIA", "TOF raw hit multiplicity - I/A side; TOF raw hits number ;Counts ",200, 0, 200) ;
155 TH1I * h2 = new TH1I("hTOFRawsOA", "TOF raw hit multiplicity - O/A side; TOF raw hits number ;Counts ",200, 0, 200) ;
156 TH1I * h3 = new TH1I("hTOFRawsIC", "TOF raw hit multiplicity - I/C side; TOF raw hits number ;Counts ",200, 0, 200) ;
157 TH1I * h4 = new TH1I("hTOFRawsOC", "TOF raw hit multiplicity - O/C side; TOF raw hits number ;Counts ",200, 0, 200) ;
159 TH1F * h5 = new TH1F("hTOFRawsTime", "TOF Raws - Hit time (ns);Measured Hit time [ns];Counts", 25000,0. ,610.) ;
161 TH1F * h6 = new TH1F("hTOFRawsTimeIA", "TOF Raws - Hit time (ns) - I/A side;Measured Hit time [ns];Counts", 25000,0. ,610.) ;
162 TH1F * h7 = new TH1F("hTOFRawsTimeOA", "TOF Raws - Hit time (ns) - O/A side;Measured Hit time [ns];Counts", 25000,0. ,610.) ;
163 TH1F * h8 = new TH1F("hTOFRawsTimeIC", "TOF Raws - Hit time (ns) - I/C side;Measured Hit time [ns];Counts", 25000,0. ,610.) ;
164 TH1F * h9 = new TH1F("hTOFRawsTimeOC", "TOF Raws - Hit time (ns) - O/C side;Measured Hit time [ns];Counts", 25000,0. ,610.) ;
166 TH1F * h10 = new TH1F("hTOFRawsToT", "TOF Raws - Hit ToT (ns);Measured Hit ToT (ns);Counts", 1000, 0., 48.8) ;
168 TH1F * h11 = new TH1F("hTOFRawsToTIA", "TOF Raws - Hit ToT (ns) - I/A side;Measured Hit ToT (ns);Counts", 1000, 0., 48.8) ;
169 TH1F * h12 = new TH1F("hTOFRawsToTOA", "TOF Raws - Hit ToT (ns) - O/A side;Measured Hit ToT (ns);Counts", 1000, 0., 48.8) ;
170 TH1F * h13 = new TH1F("hTOFRawsToTIC", "TOF Raws - Hit ToT (ns) - I/C side;Measured Hit ToT (ns);Counts", 1000, 0., 48.8) ;
171 TH1F * h14 = new TH1F("hTOFRawsToTOC", "TOF Raws - Hit ToT (ns) - O/C side;Measured Hit ToT (ns);Counts", 1000, 0., 48.8) ;
173 TH1I * h15 = new TH1I("hTOFRawsLTMHits", "LTM hits ; Crate; Counts", 72, 0., 72.);
174 TH1I * h16 = new TH1I("hTOFRawsTRMHits035", "TRM hits - crates 0 to 35 ;TRM index = SMid(crate*10)+TRM(0-9);Hits", 361, 0., 361.) ;
175 TH1I * h17 = new TH1I("hTOFRawsTRMHits3671","TRM hits - crates 36 to 71 ;TRM index = SMid(crate*10)+TRM(0-9);Hits", 361, 360., 721.) ;
177 TH1I * h18 = new TH1I("hTOFRawChannelHits","TOF channel hits count; Channel ID; Hits", 158000, 0., 158000);
179 TH1F * h19 = new TH1F("hTOFOrphansTime", "TOF Raws - Orphans time (ns);Measured Hit time [ns];Counts", 25000, 0. ,610.) ;
180 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) ;
181 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) ;
182 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) ;
183 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) ;
184 TH2F * h24 = new TH2F("hTOFTimeVsStrip","TOF Raws - Hit time vs. strip (theta); Strip index;Raws TOF time (ns) ", 91,0.,91, 250, 0., 610.) ;
186 // 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) ;
214 Add2RawsList(h0, 0, !expert, image, !saveCorr) ;
215 Add2RawsList(h1, 1, expert, !image, !saveCorr) ;
216 Add2RawsList(h2, 2, expert, !image, !saveCorr) ;
217 Add2RawsList(h3, 3, expert, !image, !saveCorr) ;
218 Add2RawsList(h4, 4, expert, !image, !saveCorr) ;
219 Add2RawsList(h5, 5, !expert, image, !saveCorr) ;
220 Add2RawsList(h6, 6, expert, !image, !saveCorr) ;
221 Add2RawsList(h7, 7, expert, !image, !saveCorr) ;
222 Add2RawsList(h8, 8, expert, !image, !saveCorr) ;
223 Add2RawsList(h9, 9, expert, !image, !saveCorr) ;
224 Add2RawsList(h10, 10, !expert, image, !saveCorr) ;
225 Add2RawsList(h11, 11, expert, !image, !saveCorr) ;
226 Add2RawsList(h12, 12, expert, !image, !saveCorr) ;
227 Add2RawsList(h13, 13, expert, !image, !saveCorr) ;
228 Add2RawsList(h14, 14, expert, !image, !saveCorr) ;
229 Add2RawsList(h15, 15, !expert, image, !saveCorr) ;
230 Add2RawsList(h16, 16, !expert, image, !saveCorr) ;
231 Add2RawsList(h17, 17, !expert, image, !saveCorr) ;
232 Add2RawsList(h18, 18, expert, !image, !saveCorr) ;
233 Add2RawsList(h19, 19, expert, !image, !saveCorr) ;
234 Add2RawsList(h20, 20, expert, !image, !saveCorr) ;
235 Add2RawsList(h21, 21, expert, !image, !saveCorr) ;
236 Add2RawsList(h22, 22, expert, !image, !saveCorr) ;
237 Add2RawsList(h23, 23, expert, !image, !saveCorr) ;
238 Add2RawsList(h24, 24, expert, !image, !saveCorr) ;
242 //____________________________________________________________________________
243 void AliTOFQADataMakerRec::InitRecPoints()
246 // create RecPoints histograms in RecPoints subdir
249 const Bool_t expert = kTRUE ;
250 const Bool_t image = kTRUE ;
252 TH1I * h0 = new TH1I("hTOFRecPoints", "TOF RecPoints multiplicity ; TOF RecPoints number;Counts",200, 0, 200) ;
254 TH1F * h1 = new TH1F("hTOFRecPointsTimeIA", "RecPoints Time Spectrum in TOF (ns)- I/A side; Calibrated TOF time [ns];Counts", 25000, 0., 610.) ;
255 TH1F * h2 = new TH1F("hTOFRecPointsTimeOA", "RecPoints Time Spectrum in TOF (ns)- O/A side; Calibrated TOF time [ns];Counts", 25000, 0., 610.) ;
256 TH1F * h3 = new TH1F("hTOFRecPointsTimeIC", "RecPoints Time Spectrum in TOF (ns)- I/C side; Calibrated TOF time [ns];Counts", 25000, 0., 610.) ;
257 TH1F * h4 = new TH1F("hTOFRecPointsTimeOC", "RecPoints Time Spectrum in TOF (ns)- O/C side; Calibrated TOF time [ns];Counts", 25000, 0., 610.) ;
259 TH1F * h5 = new TH1F("hTOFRecPointsRawTimeIA", "RecPoints raw Time Spectrum in TOF (ns)-I/A side; Measured TOF time [ns];Counts", 25000, 0., 610.) ;
260 TH1F * h6 = new TH1F("hTOFRecPointsRawTimeOA", "RecPoints raw Time Spectrum in TOF (ns)-O/A side; Measured TOF time [ns];Counts", 25000, 0., 610.) ;
261 TH1F * h7 = new TH1F("hTOFRecPointsRawTimeIC", "RecPoints raw Time Spectrum in TOF (ns)-I/C side; Measured TOF time [ns];Counts", 25000, 0., 610.) ;
262 TH1F * h8 = new TH1F("hTOFRecPointsRawTimeOC", "RecPoints raw Time Spectrum in TOF (ns)-O/C side; Measured TOF time [ns];Counts", 25000, 0., 610.) ;
264 TH1F * h9 = new TH1F("hTOFRecPointsToTIA", "RecPoints ToT Spectrum in TOF (ns)-I/A side; Measured TOT [ns];Counts", 1000, 0., 48.8 ) ;
265 TH1F * h10 = new TH1F("hTOFRecPointsToTOA", "RecPoints ToT Spectrum in TOF (ns)-O/A side; Measured TOT [ns];Counts", 1000, 0., 48.8 ) ;
266 TH1F * h11 = new TH1F("hTOFRecPointsToTIC", "RecPoints ToT Spectrum in TOF (ns)-I/C side; Measured TOT [ns];Counts", 1000, 0., 48.8 ) ;
267 TH1F * h12 = new TH1F("hTOFRecPointsToTOC", "RecPoints ToT Spectrum in TOF (ns)-O/C side; Measured TOT [ns];Counts", 1000, 0., 48.8 ) ;
268 TH1F * h13 = new TH1F("hTOFRawChannelHits","TOF channel hits count; Channel ID; Hits", 158000, 0., 158000);
270 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) ;
271 TH2F * h15 = new TH2F("hTOFRecTimeVsStrip","RecPoints TOF time vs. strip (theta); Strip index; RecPoints TOF time (ns) ",92,-1.,91, 250, 0., 610.) ;
273 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) ;
274 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) ;
275 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) ;
276 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) ;
299 Add2RecPointsList(h0, 0, !expert, image) ;
300 Add2RecPointsList(h1, 1, !expert, image) ;
301 Add2RecPointsList(h2, 2, !expert, image) ;
302 Add2RecPointsList(h3, 3, !expert, image) ;
303 Add2RecPointsList(h4, 4, !expert, image) ;
304 Add2RecPointsList(h5, 5, expert, !image) ;
305 Add2RecPointsList(h6, 6, expert, !image) ;
306 Add2RecPointsList(h7, 7, expert, !image) ;
307 Add2RecPointsList(h8, 8, expert, !image) ;
308 Add2RecPointsList(h9, 9, !expert, image) ;
309 Add2RecPointsList(h10, 10, !expert, image) ;
310 Add2RecPointsList(h11, 11, !expert, image) ;
311 Add2RecPointsList(h12, 12, !expert, image) ;
312 Add2RecPointsList(h13, 13, expert, !image) ;
313 Add2RecPointsList(h14, 14, expert, !image) ;
314 Add2RecPointsList(h15, 15, expert, !image) ;
315 Add2RecPointsList(h16, 16, expert, !image) ;
316 Add2RecPointsList(h17, 17, expert, !image) ;
317 Add2RecPointsList(h18, 18, expert, !image) ;
318 Add2RecPointsList(h19, 19, expert, !image) ;
321 //____________________________________________________________________________
322 void AliTOFQADataMakerRec::InitESDs()
325 //create ESDs histograms in ESDs subdir
328 const Bool_t expert = kTRUE ;
329 const Bool_t image = kTRUE ;
331 TH1I * h0 = new TH1I("hTOFESDs", "Number of matched TOF tracks per event;Number of TOF matched ESD tracks;Counts", 200, 0, 200) ;
332 TH1F * h1 = new TH1F("hTOFESDsTime", "Matched ESDs tracks: TOF Time spectrum; Calibrated TOF time [ns];Counts", 25000, 0., 610. ) ;
333 TH1F * h2 = new TH1F("hTOFESDsRawTime", "Matched ESDs tracks: TOF raw Time spectrum;Measured TOF time [ns];Counts", 25000, 0., 610.) ;
334 TH1F * h3 = new TH1F("hTOFESDsToT", "Matched ESDs tracks: TOF ToT spectrum; ESDs ToT [ns];Counts",1000, 0., 48.8) ;
335 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.) ;
336 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.) ;
337 TH1F * h6 = new TH1F("hTOFESDsTPCmatched", "TPC-TOF matched tracks momentum distribution (GeV/c); p (GeV/c);Counts", 50, 0.20, 5.00) ;
338 TH1F * h7 = new TH1F("hTOFESDsMatchingProb", "TPC-TOF track-matching probability per event;TPC-TOF track-matching probability (%) ;Counts",101, -1.0, 100) ;
339 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) ;
340 TH1F * h9 = new TH1F("hTOFESDsDiffTime", "ESDs t_{TOF}-t_{exp} spectrum in TOF (ps); t_{TOF}-t_{exp} [ps];Counts", 200, -2440., 2440.) ;
341 TH1F * h10 = new TH1F("hTOFHitsLength", "Matched ESDs tracks: Length Spectrum; Track length [cm];Counts", 800, 0., 800) ;
354 Add2ESDsList(h0, 0, !expert, image) ;
355 Add2ESDsList(h1, 1, !expert, image) ;
356 Add2ESDsList(h2, 2, expert, image) ;
357 Add2ESDsList(h3, 3, !expert, image) ;
358 Add2ESDsList(h4, 4, expert, image) ;
359 Add2ESDsList(h5, 5, expert, image) ;
360 Add2ESDsList(h6, 6, expert, image) ;
361 Add2ESDsList(h7, 7, expert, image) ;
362 Add2ESDsList(h8, 8, expert, image) ;
363 Add2ESDsList(h9, 9, !expert, image) ;
364 Add2ESDsList(h10, 10, !expert, image) ;
368 //____________________________________________________________________________
369 void AliTOFQADataMakerRec::MakeRaws(AliRawReader* rawReader)
372 // makes data from Raws
374 if (rawReader->GetType()==7) {
376 Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;//in ns
377 Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
379 Int_t ntof[5]; /* 0=tot, 1=IA, 2=OA, 3=IC, 4=OC*/
380 for (Int_t j=0;j<5;j++){ ntof[j]=0;}
382 Int_t equipmentID[5]; //(ddl, trm, chain,tdc,channel)
383 Int_t volumeID[5]; //(sector,plate,strip,padX,padZ)
384 Int_t volumeID2[5]; //(sector,plate,strip,padZ,padX) to use AliTOFGeometry::GetIndex()
385 Int_t out[5]; // out=(indexZ,indexPhi)
388 TClonesArray * clonesRawData;
389 AliTOFRawStream tofInput(rawReader);
391 //uncomment if needed to apply DeltaBC correction
392 //tofInput.ApplyBCCorrections(kTRUE);
394 for (Int_t iDDL = 0; iDDL < AliTOFGeometry::NDDL()*AliTOFGeometry::NSectors(); iDDL++){
397 tofInput.LoadRawDataBuffers(iDDL);
398 clonesRawData = (TClonesArray*)tofInput.GetRawData();
399 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
400 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
402 if (tofRawDatum->GetTOF()){
405 equipmentID[1]=tofRawDatum->GetTRM();
406 equipmentID[2]=tofRawDatum->GetTRMchain();
407 equipmentID[3]=tofRawDatum->GetTDC();
408 equipmentID[4]=tofRawDatum->GetTDCchannel();
410 if (CheckEquipID(equipmentID)){
411 tofInput.EquipmentId2VolumeId(iDDL,
412 tofRawDatum->GetTRM(),
413 tofRawDatum->GetTRMchain(),
414 tofRawDatum->GetTDC(),
415 tofRawDatum->GetTDCchannel(),
417 if (FilterSpare(equipmentID)) continue;
418 if (FilterLTMData(equipmentID)){ //counts LTM hits
419 if (tofRawDatum->GetTOT()) GetRawsData(15)->Fill(equipmentID[0]);
421 if (CheckVolumeID(volumeID)){
423 GetMapIndeces(volumeID,out);
424 volumeID2[0]=volumeID[0];
425 volumeID2[1]=volumeID[1];
426 volumeID2[2]=volumeID[2];
427 volumeID2[3]=volumeID[4];
428 volumeID2[4]=volumeID[3];
429 chIndex=AliTOFGeometry::GetIndex(volumeID2);
431 if (tofRawDatum->GetTOT()){
432 if (!(fCalibData->GetNoiseStatus(chIndex)==AliTOFChannelOnlineStatusArray::kTOFNoiseBad)
433 && (fCalibData->GetHWStatus(chIndex) == AliTOFChannelOnlineStatusArray::kTOFHWOk)) {//noise and enabled filter
434 ntof[0]++; //counter for tof hits
436 //fill global spectra for DQM plots
437 GetRawsData(5)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;//in ns
438 GetRawsData(10)->Fill( tofRawDatum->GetTOT()*tot2ns) ;//in ns
440 //fill side-related spectra for experts plots
441 if (volumeID2[0]>4 && volumeID2[0]<14){ //I side
442 if ((iDDL%4==0)|| (iDDL%4==1)){ //A side
444 GetRawsData(6)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;
445 GetRawsData(11)->Fill( tofRawDatum->GetTOT()*tot2ns) ;
447 if ((iDDL%4==2)|| (iDDL%4==3)){//C side
449 GetRawsData(8)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;
450 GetRawsData(13)->Fill( tofRawDatum->GetTOT()*tot2ns) ;
454 if (volumeID2[0]<5 || volumeID2[0]>13){ //O side
455 if ((iDDL%4==0)|| (iDDL%4==1)){ //A side
457 GetRawsData(7)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;
458 GetRawsData(12)->Fill( tofRawDatum->GetTOT()*tot2ns) ;
460 if ((iDDL%4==2)|| (iDDL%4==3)){//C side
462 GetRawsData(9)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;
463 GetRawsData(14)->Fill( tofRawDatum->GetTOT()*tot2ns) ;
468 GetRawsData(18)->Fill(chIndex);
470 Int_t trm= iDDL*10+(equipmentID[1]-3);
471 if (iDDL>=0 && iDDL<36) {
472 GetRawsData(16)->Fill(trm) ;
473 GetRawsData(20)->Fill(trm,tofRawDatum->GetTOF()*tdc2ns);
474 GetRawsData(22)->Fill(trm,tofRawDatum->GetTOT()*tot2ns);
476 if (iDDL>=36 && iDDL<72) {
477 GetRawsData(17)->Fill(trm) ;//in ns
478 GetRawsData(21)->Fill(trm,tofRawDatum->GetTOF()*tdc2ns);
479 GetRawsData(23)->Fill(trm,tofRawDatum->GetTOT()*tot2ns);
481 GetRawsData(24)->Fill(GetStripIndex(volumeID),tofRawDatum->GetTOF()*tdc2ns) ;
482 //GetRawsData(25)->Fill( out[0],out[1]) ;//raw map
486 if (!(fCalibData->GetNoiseStatus(chIndex) == AliTOFChannelOnlineStatusArray::kTOFNoiseBad)
487 && (fCalibData->GetHWStatus(chIndex) == AliTOFChannelOnlineStatusArray::kTOFHWOk))
488 GetRawsData(19)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;//in ns
490 }//end volumeID check
495 clonesRawData->Clear();
498 for (Int_t j=0;j<5;j++){
499 GetRawsData(j)->Fill(ntof[j]);
501 fProcessedRawEventN++;
504 AliDebug(1,Form("Event of type %d found. Skipping non-physics event for QA.\n", rawReader->GetType()));
506 EnableDqmShifterOpt(kTRUE);
509 //____________________________________________________________________________
510 void AliTOFQADataMakerRec::MakeRecPoints(TTree * clustersTree)
513 // Make data from Clusters
516 Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
517 Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
523 TBranch *branch=clustersTree->GetBranch("TOF");
525 AliError("can't get the branch with the TOF clusters !");
529 static TClonesArray dummy("AliTOFcluster",10000);
531 TClonesArray *clusters=&dummy;
532 branch->SetAddress(&clusters);
535 clustersTree->GetEvent(0);
537 GetRecPointsData(0)->Fill((Int_t)clusters->GetEntriesFast()) ;
539 TIter next(clusters) ;
541 while ( (c = dynamic_cast<AliTOFcluster *>(next())) ) {
543 volumeID[0] = c->GetDetInd(0);
544 volumeID[1] = c->GetDetInd(1);
545 volumeID[2] = c->GetDetInd(2);
546 volumeID[3] = c->GetDetInd(4); //X and Z indeces inverted in RecPoints
547 volumeID[4] = c->GetDetInd(3); //X and Z indeces inverted in RecPoints
549 for (Int_t i=0;i<5;i++){
550 volumeID2[i]=c->GetDetInd(i); //X and Z indeces inverted in RecPoints
553 GetMapIndeces(volumeID,out);
554 Int_t chIndex=AliTOFGeometry::GetIndex(volumeID2);
555 Int_t iDDL=AliTOFRawStream::Geant2DDL(volumeID2);
556 Int_t iTRM=AliTOFRawStream::Geant2TRM(volumeID2);
558 if (c->GetTDCRAW() && c->GetTDC() && c->GetToT()){
559 if (volumeID2[0]>4 && volumeID2[0]<14){ //I side
560 if ((iDDL%4==0)|| (iDDL%4==1)){ //A side
561 GetRecPointsData(1)->Fill( c->GetTDC()*tdc2ns) ;//in ns
562 GetRecPointsData(5)->Fill( c->GetTDCRAW()*tdc2ns) ;//in ns
563 GetRecPointsData(9)->Fill( c->GetToT()*tot2ns) ;//in ns
566 if ((iDDL%4==2)|| (iDDL%4==3)){//C side
567 GetRecPointsData(3)->Fill( c->GetTDC()*tdc2ns) ;//in ns
568 GetRecPointsData(7)->Fill( c->GetTDCRAW()*tdc2ns) ;//in ns
569 GetRecPointsData(11)->Fill( c->GetToT()*tot2ns) ;//in ns
573 if (volumeID2[0]<5 || volumeID2[0]>13){ //O side
574 if ((iDDL%4==0)|| (iDDL%4==1)){ //A side
575 GetRecPointsData(2)->Fill( c->GetTDC()*tdc2ns) ;//in ns
576 GetRecPointsData(6)->Fill( c->GetTDCRAW()*tdc2ns) ;//in ns
577 GetRecPointsData(10)->Fill( c->GetToT()*tot2ns) ;//in ns
579 if ((iDDL%4==2)|| (iDDL%4==3)){//C side
580 GetRecPointsData(4)->Fill( c->GetTDC()*tdc2ns) ;//in ns
581 GetRecPointsData(8)->Fill( c->GetTDCRAW()*tdc2ns) ;//in ns
582 GetRecPointsData(12)->Fill( c->GetToT()*tot2ns) ;//in ns
587 GetRecPointsData(14)->Fill(out[0],out[1]);
588 GetRecPointsData(15)->Fill(GetStripIndex(volumeID), c->GetTDC()*tdc2ns) ;
589 Int_t trm= iDDL*10+(iTRM-3);
590 if (iDDL>=0 && iDDL<36) {
591 GetRecPointsData(16)->Fill(trm,c->GetTDC()*tdc2ns);
592 GetRecPointsData(18)->Fill(trm,c->GetToT()*tot2ns);
594 if (iDDL>=36 && iDDL<72) {
595 GetRecPointsData(17)->Fill(trm,c->GetTDC()*tdc2ns);
596 GetRecPointsData(19)->Fill(trm,c->GetToT()*tot2ns);
598 GetRecPointsData(13)->Fill(chIndex) ;//in ns
602 EnableDqmShifterOpt(kFALSE);
605 //____________________________________________________________________________
606 void AliTOFQADataMakerRec::MakeESDs(AliESDEvent * esd)
609 // make QA data from ESDs
611 const Double_t speedOfLight = TMath::C()*1E2*1E-12; // cm/ps
612 const Double_t pionMass = 0.13957018; //GeV/c^2
614 Int_t ntrk = esd->GetNumberOfTracks() ;
623 AliESDtrack *track=esd->GetTrack(ntrk);
624 Double_t tofTime=track->GetTOFsignal();//in ps
625 Double_t tofTimeRaw=track->GetTOFsignalRaw();//in ps
626 Double_t tofToT=track->GetTOFsignalToT(); //in ps
628 UInt_t status=track->GetStatus();
629 if (track->IsOn(AliESDtrack::kTPCrefit)) {
631 Double_t y=track->Eta();
632 if ((status&AliESDtrack::kTOFout)!=0) {
637 if (TMath::Abs(y)<0.9) {
638 GetESDsData(6)->Fill(track->Pt());
641 GetESDsData(1)->Fill(tofTime*1E-3);
642 GetESDsData(2)->Fill(tofTimeRaw*1E-3);
643 GetESDsData(3)->Fill(tofToT*1E-3);
644 //check how many tracks where ESD PID is ok
645 if ((status&AliESDtrack::kESDpid)!=0) nesdpid++;
646 if (((status&AliESDtrack::kESDpid)&AliESDtrack::kTOFpid)!=0) ntofpid++;
648 Double_t length =track->GetIntegratedLength();
649 Double_t mom2=(track->Pt()*track->Pt())+(track->Pz()*track->Pz());
650 Double_t eTexp = TMath::Sqrt(1+(pionMass*pionMass/mom2))*length/speedOfLight; //in ps
651 GetESDsData(9)->Fill(tofTime-eTexp);
652 GetESDsData(10)->Fill(length);
653 } //end check on matched tracks
655 }//end check on TPCrefit
658 GetESDsData(0)->Fill(ntof) ;
662 Float_t ratio = (Int_t)ntofpid/(Int_t)ntof*100; //identified by TOF over matched
663 GetESDsData(4)->Fill(ratio) ;
667 Float_t ratio = (Float_t)ntofpid/(Float_t)nesdpid *100; //identified by TOF over identified
668 GetESDsData(5)->Fill(ratio) ;
672 Float_t ratio = (Float_t)ntof/(Float_t)ntpc*100.; //matching probability
673 GetESDsData(7)->Fill(ratio) ;
677 Float_t ratio = (Float_t)ntof/(Float_t)ntofout*100; //matched over propagated to TOF outer radius
678 GetESDsData(8)->Fill(ratio) ;
680 EnableDqmShifterOpt(kFALSE);
683 //____________________________________________________________________________
684 void AliTOFQADataMakerRec::StartOfDetectorCycle()
687 //Detector specific actions at start of cycle
689 fCalibData = GetCalibData();
693 //____________________________________________________________________________
694 void AliTOFQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
696 //Detector specific actions at end of cycle
697 // do the QA checking
698 for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
699 if ( !AliQAv1::Instance()->IsEventSpecieSet(specie) )
702 AliInfo(Form("Processed %i physics raw events",fProcessedRawEventN));
704 if (fEnableDqmShifterOpt){
705 // Help make the raw qa histogram easier to interpret for the DQM shifter
706 // This is still to be optimized...
708 if (!GetRawsData(0) || !GetRawsData(5) || !GetRawsData(10)
709 || !GetRawsData(15) || !GetRawsData(16) || !GetRawsData(17)) {
710 printf("No histogram for DQM found - Possible memory corruption ???. Please check\n") ;
714 Double_t yTimeMin = GetRawsData(5)->GetMinimum();
715 Double_t yTimeMax = GetRawsData(5)->GetMaximum();
716 Double_t yTotMin = GetRawsData(10)->GetMinimum();
717 Double_t yTotMax = GetRawsData(10)->GetMaximum();
719 TLine* lineExpTimeMin = new TLine(200., yTimeMin, 200., yTimeMax);
720 lineExpTimeMin->SetLineColor(kGreen);
721 lineExpTimeMin->SetLineWidth(2);
723 TLine* lineExpTimeMax = new TLine(250., yTimeMin, 250., yTimeMax);
724 lineExpTimeMax->SetLineColor(kGreen);
725 lineExpTimeMax->SetLineWidth(2);
727 TLine* lineExpTotMin = new TLine( 5., yTotMin, 5., yTotMax);
728 lineExpTotMin->SetLineColor(kGreen);
729 lineExpTotMin->SetLineWidth(2);
731 TLine* lineExpTotMax = new TLine(20., yTotMin, 20., yTotMax);
732 lineExpTotMax->SetLineColor(kGreen);
733 lineExpTotMax->SetLineWidth(2);
735 ((TH1F*)GetRawsData(5))->GetListOfFunctions()->Add(lineExpTimeMin);
736 ((TH1F*)GetRawsData(5))->GetListOfFunctions()->Add(lineExpTimeMax);
737 ((TH1F*)GetRawsData(10))->GetListOfFunctions()->Add(lineExpTotMin);
738 ((TH1F*)GetRawsData(10))->GetListOfFunctions()->Add(lineExpTotMax);
740 //make up for all bistos
741 for(Int_t j=0;j<20;j++){
743 GetRawsData(j)->SetMarkerColor(kRed);
744 GetRawsData(j)->SetMarkerStyle(7);
746 GetRawsData(j)->SetLineColor(kBlue+1);
747 GetRawsData(j)->SetLineWidth(1);
748 GetRawsData(j)->SetMarkerColor(kBlue+1);
751 Float_t ySMmax035=GetRawsData(16)->GetMaximum();
752 TLine* lineSMid035[10];
753 Float_t ySMmax3671=GetRawsData(17)->GetMaximum();
754 TLine* lineSMid3671[10];
756 for (Int_t sm=0;sm<10;sm++){
757 lineSMid035[sm] = new TLine( 40*sm, 0, 40*sm, ySMmax035);
758 lineSMid035[sm]->SetLineColor(kMagenta);
759 lineSMid035[sm]->SetLineWidth(2);
760 GetRawsData(16)->GetListOfFunctions()->Add(lineSMid035[sm]);
762 lineSMid3671[sm] = new TLine( 40*sm+360, 0, 40*sm+360, ySMmax3671);
763 lineSMid3671[sm]->SetLineColor(kMagenta);
764 lineSMid3671[sm]->SetLineWidth(2);
765 GetRawsData(17)->GetListOfFunctions()->Add(lineSMid3671[sm]);
768 for (Int_t j=15;j<19;j++){
769 GetRawsData(j)->SetFillColor(kGray+1);
770 GetRawsData(j)->SetOption("bar");
775 AliQAChecker::Instance()->Run(AliQAv1::kTOF, task, list) ;
777 //____________________________________________________________________________
778 void AliTOFQADataMakerRec::GetMapIndeces(const Int_t* const in , Int_t* out)
781 //return appropriate indeces for the theta-phi map
784 Int_t npadX = AliTOFGeometry::NpadX();
785 Int_t npadZ = AliTOFGeometry::NpadZ();
786 Int_t nStripA = AliTOFGeometry::NStripA();
787 Int_t nStripB = AliTOFGeometry::NStripB();
788 Int_t nStripC = AliTOFGeometry::NStripC();
790 Int_t isector = in[0];
791 Int_t iplate = in[1];
792 Int_t istrip = in[2];
796 Int_t stripOffset = 0;
802 stripOffset = nStripC;
805 stripOffset = nStripC+nStripB;
808 stripOffset = nStripC+nStripB+nStripA;
811 stripOffset = nStripC+nStripB+nStripA+nStripB;
814 AliDebug(1,Form("Wrong plate number in TOF (%d) !",iplate));
817 Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
818 Int_t phiindex=npadX*isector+ipadX+1;
824 //---------------------------------------------------------------
825 Int_t AliTOFQADataMakerRec::GetStripIndex(const Int_t * const in)
827 /* return tof strip index between 0 and 91 */
829 Int_t nStripA = AliTOFGeometry::NStripA();
830 Int_t nStripB = AliTOFGeometry::NStripB();
831 Int_t nStripC = AliTOFGeometry::NStripC();
833 // Int_t isector = in[0];
834 Int_t iplate = in[1];
835 Int_t istrip = in[2];
836 //Int_t ipadX = in[3];
837 //Int_t ipadZ = in[4];
839 Int_t stripOffset = 0;
845 stripOffset = nStripC;
848 stripOffset = nStripC+nStripB;
851 stripOffset = nStripC+nStripB+nStripA;
854 stripOffset = nStripC+nStripB+nStripA+nStripB;
857 AliDebug(1,Form("Wrong plate number in TOF (%d) !",iplate));
862 if (stripOffset<0 || stripOffset>92) return -1;
864 return (stripOffset+istrip);
867 //---------------------------------------------------------------
868 Bool_t AliTOFQADataMakerRec::CheckVolumeID(const Int_t * const volumeID)
871 //Checks volume ID validity
874 for (Int_t j=0;j<5;j++){
876 AliDebug(1,Form("Invalid detector volume index for volumeID[%i]",j));
884 //---------------------------------------------------------------
885 Bool_t AliTOFQADataMakerRec::CheckEquipID(const Int_t * const equipmentID)
888 //Checks equipment ID validity
890 for (Int_t j=0;j<5;j++){
891 if (equipmentID[j]<0) {
892 AliDebug(1,Form("Invalid equipment volume index for equipmentID[%i]",j));
898 //---------------------------------------------------------------
899 Bool_t AliTOFQADataMakerRec::FilterLTMData(const Int_t * const equipmentID) const
901 /*It returns kTRUE if data come from LTM.
902 It thus filters trigger-related signals */
905 //if (!CheckEquipID(equipmentID)) return kFALSE;
906 ddl = equipmentID[0];
907 trm = equipmentID[1];
908 tdc = equipmentID[3];
910 if ((ddl%2==1) && (trm==3) && (tdc>11 && tdc<15))
916 //---------------------------------------------------------------
917 Bool_t AliTOFQADataMakerRec::FilterSpare(const Int_t * const equipmentID) const
919 /*It returns kTRUE if data come from spare
921 So far only check on TRM 3 crate left is implemented */
924 //if (!CheckEquipID(equipmentID)) return kFALSE;
925 ddl = equipmentID[0];
926 trm = equipmentID[1];
927 tdc = equipmentID[3];
929 if ((ddl%2==1) && (trm==3) && (tdc>2 && tdc<12))