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 ///////////////////////////////////////////////////////////////////////
24 /* Modified by fbellini on 22/04/2010
25 - Added filter for physics events
27 Modified by fbellini on 16/04/2010
28 - Added EnableDqmShifterOpt()
29 - Modified EndOfDetectorCycle() with options for DQM
32 Modified by fbellini on 30/03/2010
33 - Changed raws time histos range to 610ns
34 - Added FilterLTMData() and FilterSpare() methods
35 - Added check on enabled channels for raw data
36 - Updated RecPoints QA
38 Modified by fbellini on 02/03/2010
39 - Fixed raw data decoding methods (use AliTOFRawStream::LoadRawDataBuffer())
40 - Added filter for noisy channels and read map from OCDB
41 - Added GetCalibData() method
42 - Added CheckVolumeID() and CheckEquipID() methods
46 #include <TClonesArray.h>
51 #include "AliCDBManager.h"
52 #include "AliCDBEntry.h"
53 #include "AliESDEvent.h"
54 #include "AliESDtrack.h"
55 #include "AliQAChecker.h"
56 #include "AliRawReader.h"
57 #include "AliTOFRawStream.h"
58 #include "AliTOFcluster.h"
59 #include "AliTOFQADataMakerRec.h"
60 #include "AliTOFrawData.h"
61 #include "AliTOFGeometry.h"
62 #include "AliTOFChannelOnlineStatusArray.h"
65 ClassImp(AliTOFQADataMakerRec)
67 //____________________________________________________________________________
68 AliTOFQADataMakerRec::AliTOFQADataMakerRec() :
69 AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kTOF), "TOF Quality Assurance Data Maker"),
71 fEnableNoiseFiltering(kFALSE),
72 fEnableDqmShifterOpt(kFALSE),
73 fProcessedRawEventN(0)
81 //____________________________________________________________________________
82 AliTOFQADataMakerRec::AliTOFQADataMakerRec(const AliTOFQADataMakerRec& qadm) :
84 fCalibData(qadm.fCalibData),
85 fEnableNoiseFiltering(qadm.fEnableNoiseFiltering),
86 fEnableDqmShifterOpt(qadm.fEnableDqmShifterOpt),
87 fProcessedRawEventN(qadm.fProcessedRawEventN)
92 SetName((const char*)qadm.GetName()) ;
93 SetTitle((const char*)qadm.GetTitle());
98 //__________________________________________________________________
99 AliTOFQADataMakerRec& AliTOFQADataMakerRec::operator = (const AliTOFQADataMakerRec& qadm )
102 // assignment operator.
104 this->~AliTOFQADataMakerRec();
105 new(this) AliTOFQADataMakerRec(qadm);
109 //----------------------------------------------------------------------------
110 AliTOFChannelOnlineStatusArray* AliTOFQADataMakerRec::GetCalibData() const
113 // Retrive TOF calib objects from OCDB
115 AliCDBManager *man = AliCDBManager::Instance();
119 cdbe = man->Get("TOF/Calib/Status",fRun);
121 AliWarning("Load of calibration data from default storage failed!");
122 AliWarning("Calibration data will be loaded from local storage ($ALICE_ROOT)");
123 man->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
124 cdbe = man->Get("TOF/Calib/Status",fRun);
126 // Retrieval of data in directory TOF/Calib/Data:
128 AliTOFChannelOnlineStatusArray * array = 0;
129 if (cdbe) array = (AliTOFChannelOnlineStatusArray *)cdbe->GetObject();
130 if (!array) AliFatal("No calibration data from calibration database !");
138 //____________________________________________________________________________
139 void AliTOFQADataMakerRec::InitRaws()
142 // create Raws histograms in Raws subdir
144 const Bool_t expert = kTRUE ;
145 const Bool_t saveCorr = kTRUE ;
146 const Bool_t image = kTRUE ;
148 TH1F * h0 = new TH1F("hTOFRaws", "TOF raw hit multiplicity; TOF raw hits number;Counts ",201, -1., 200.) ;
150 TH1F * h1 = new TH1F("hTOFRawsIA", "TOF raw hit multiplicity - I/A side; TOF raw hits number ;Counts ",201, -1., 200.) ;
151 TH1F * h2 = new TH1F("hTOFRawsOA", "TOF raw hit multiplicity - O/A side; TOF raw hits number ;Counts ",201, -1., 200.) ;
152 TH1F * h3 = new TH1F("hTOFRawsIC", "TOF raw hit multiplicity - I/C side; TOF raw hits number ;Counts ",201, -1., 200.) ;
153 TH1F * h4 = new TH1F("hTOFRawsOC", "TOF raw hit multiplicity - O/C side; TOF raw hits number ;Counts ",201, -1., 200.) ;
155 TH1F * h5 = new TH1F("hTOFRawsTime", "TOF Raws - Hit time (ns);Measured Hit time [ns];Counts", 25000,0. ,610.) ;
157 TH1F * h6 = new TH1F("hTOFRawsTimeIA", "TOF Raws - Hit time (ns) - I/A side;Measured Hit time [ns];Counts", 25000,0. ,610.) ;
158 TH1F * h7 = new TH1F("hTOFRawsTimeOA", "TOF Raws - Hit time (ns) - O/A side;Measured Hit time [ns];Counts", 25000,0. ,610.) ;
159 TH1F * h8 = new TH1F("hTOFRawsTimeIC", "TOF Raws - Hit time (ns) - I/C side;Measured Hit time [ns];Counts", 25000,0. ,610.) ;
160 TH1F * h9 = new TH1F("hTOFRawsTimeOC", "TOF Raws - Hit time (ns) - O/C side;Measured Hit time [ns];Counts", 25000,0. ,610.) ;
162 TH1F * h10 = new TH1F("hTOFRawsToT", "TOF Raws - Hit ToT (ns);Measured Hit ToT (ns);Counts", 1000, 0., 48.8) ;
164 TH1F * h11 = new TH1F("hTOFRawsToTIA", "TOF Raws - Hit ToT (ns) - I/A side;Measured Hit ToT (ns);Counts", 1000, 0., 48.8) ;
165 TH1F * h12 = new TH1F("hTOFRawsToTOA", "TOF Raws - Hit ToT (ns) - O/A side;Measured Hit ToT (ns);Counts", 1000, 0., 48.8) ;
166 TH1F * h13 = new TH1F("hTOFRawsToTIC", "TOF Raws - Hit ToT (ns) - I/C side;Measured Hit ToT (ns);Counts", 1000, 0., 48.8) ;
167 TH1F * h14 = new TH1F("hTOFRawsToTOC", "TOF Raws - Hit ToT (ns) - O/C side;Measured Hit ToT (ns);Counts", 1000, 0., 48.8) ;
169 TH1F * h15 = new TH1F("hTOFRawsLTMHits", "LTM hits ; Crate; Counts", 72, 0., 72.);
170 TH1F * h16 = new TH1F("hTOFRawsTRMHits035", "TRM hits - crates 0 to 35 ;TRM index = SMid(crate*10)+TRM(0-9);Hits", 361, 0., 361.) ;
171 TH1F * h17 = new TH1F("hTOFRawsTRMHits3671","TRM hits - crates 36 to 71 ;TRM index = SMid(crate*10)+TRM(0-9);Hits", 361, 360., 721.) ;
173 TH1F * h18 = new TH1F("hTOFRawChannelHits","TOF channel hits count; Channel ID; Hits", 158000, 0., 158000);
175 TH1F * h19 = new TH1F("hTOFOrphansTime", "TOF Raws - Orphans time (ns);Measured Hit time [ns];Counts", 25000, 0. ,610.) ;
176 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) ;
177 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) ;
178 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) ;
179 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) ;
180 TH2F * h24 = new TH2F("hTOFTimeVsStrip","TOF Raws - Hit time vs. strip (theta); Strip index;Raws TOF time (ns) ", 91,0.,91, 250, 0., 610.) ;
182 // 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) ;
210 Add2RawsList(h0, 0, !expert, image, !saveCorr) ;
211 Add2RawsList(h1, 1, expert, !image, !saveCorr) ;
212 Add2RawsList(h2, 2, expert, !image, !saveCorr) ;
213 Add2RawsList(h3, 3, expert, !image, !saveCorr) ;
214 Add2RawsList(h4, 4, expert, !image, !saveCorr) ;
215 Add2RawsList(h5, 5, !expert, image, !saveCorr) ;
216 Add2RawsList(h6, 6, expert, !image, !saveCorr) ;
217 Add2RawsList(h7, 7, expert, !image, !saveCorr) ;
218 Add2RawsList(h8, 8, expert, !image, !saveCorr) ;
219 Add2RawsList(h9, 9, expert, !image, !saveCorr) ;
220 Add2RawsList(h10, 10, !expert, image, !saveCorr) ;
221 Add2RawsList(h11, 11, expert, !image, !saveCorr) ;
222 Add2RawsList(h12, 12, expert, !image, !saveCorr) ;
223 Add2RawsList(h13, 13, expert, !image, !saveCorr) ;
224 Add2RawsList(h14, 14, expert, !image, !saveCorr) ;
225 Add2RawsList(h15, 15, !expert, image, !saveCorr) ;
226 Add2RawsList(h16, 16, !expert, image, !saveCorr) ;
227 Add2RawsList(h17, 17, !expert, image, !saveCorr) ;
228 Add2RawsList(h18, 18, expert, !image, !saveCorr) ;
229 Add2RawsList(h19, 19, expert, !image, !saveCorr) ;
230 Add2RawsList(h20, 20, expert, !image, !saveCorr) ;
231 Add2RawsList(h21, 21, expert, !image, !saveCorr) ;
232 Add2RawsList(h22, 22, expert, !image, !saveCorr) ;
233 Add2RawsList(h23, 23, expert, !image, !saveCorr) ;
234 Add2RawsList(h24, 24, expert, !image, !saveCorr) ;
238 //____________________________________________________________________________
239 void AliTOFQADataMakerRec::InitRecPoints()
242 // create RecPoints histograms in RecPoints subdir
245 const Bool_t expert = kTRUE ;
246 const Bool_t image = kTRUE ;
248 TH1F * h0 = new TH1F("hTOFRecPoints", "TOF RecPoints multiplicity ; TOF RecPoints number;Counts",201, -1., 200.) ;
250 TH1F * h1 = new TH1F("hTOFRecPointsTimeIA", "RecPoints Time Spectrum in TOF (ns)- I/A side; Calibrated TOF time [ns];Counts", 25000, 0., 610.) ;
251 TH1F * h2 = new TH1F("hTOFRecPointsTimeOA", "RecPoints Time Spectrum in TOF (ns)- O/A side; Calibrated TOF time [ns];Counts", 25000, 0., 610.) ;
252 TH1F * h3 = new TH1F("hTOFRecPointsTimeIC", "RecPoints Time Spectrum in TOF (ns)- I/C side; Calibrated TOF time [ns];Counts", 25000, 0., 610.) ;
253 TH1F * h4 = new TH1F("hTOFRecPointsTimeOC", "RecPoints Time Spectrum in TOF (ns)- O/C side; Calibrated TOF time [ns];Counts", 25000, 0., 610.) ;
255 TH1F * h5 = new TH1F("hTOFRecPointsRawTimeIA", "RecPoints raw Time Spectrum in TOF (ns)-I/A side; Measured TOF time [ns];Counts", 25000, 0., 610.) ;
256 TH1F * h6 = new TH1F("hTOFRecPointsRawTimeOA", "RecPoints raw Time Spectrum in TOF (ns)-O/A side; Measured TOF time [ns];Counts", 25000, 0., 610.) ;
257 TH1F * h7 = new TH1F("hTOFRecPointsRawTimeIC", "RecPoints raw Time Spectrum in TOF (ns)-I/C side; Measured TOF time [ns];Counts", 25000, 0., 610.) ;
258 TH1F * h8 = new TH1F("hTOFRecPointsRawTimeOC", "RecPoints raw Time Spectrum in TOF (ns)-O/C side; Measured TOF time [ns];Counts", 25000, 0., 610.) ;
260 TH1F * h9 = new TH1F("hTOFRecPointsToTIA", "RecPoints ToT Spectrum in TOF (ns)-I/A side; Measured TOT [ns];Counts", 1000, 0., 48.8 ) ;
261 TH1F * h10 = new TH1F("hTOFRecPointsToTOA", "RecPoints ToT Spectrum in TOF (ns)-O/A side; Measured TOT [ns];Counts", 1000, 0., 48.8 ) ;
262 TH1F * h11 = new TH1F("hTOFRecPointsToTIC", "RecPoints ToT Spectrum in TOF (ns)-I/C side; Measured TOT [ns];Counts", 1000, 0., 48.8 ) ;
263 TH1F * h12 = new TH1F("hTOFRecPointsToTOC", "RecPoints ToT Spectrum in TOF (ns)-O/C side; Measured TOT [ns];Counts", 1000, 0., 48.8 ) ;
264 TH1F * h13 = new TH1F("hTOFRawChannelHits","TOF channel hits count; Channel ID; Hits", 158000, 0., 158000);
266 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) ;
267 TH2F * h15 = new TH2F("hTOFRecTimeVsStrip","RecPoints TOF time vs. strip (theta); Strip index; RecPoints TOF time (ns) ",92,-1.,91, 250, 0., 610.) ;
269 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) ;
270 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) ;
271 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) ;
272 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) ;
295 Add2RecPointsList(h0, 0, !expert, image) ;
296 Add2RecPointsList(h1, 1, !expert, image) ;
297 Add2RecPointsList(h2, 2, !expert, image) ;
298 Add2RecPointsList(h3, 3, !expert, image) ;
299 Add2RecPointsList(h4, 4, !expert, image) ;
300 Add2RecPointsList(h5, 5, expert, !image) ;
301 Add2RecPointsList(h6, 6, expert, !image) ;
302 Add2RecPointsList(h7, 7, expert, !image) ;
303 Add2RecPointsList(h8, 8, expert, !image) ;
304 Add2RecPointsList(h9, 9, !expert, image) ;
305 Add2RecPointsList(h10, 10, !expert, image) ;
306 Add2RecPointsList(h11, 11, !expert, image) ;
307 Add2RecPointsList(h12, 12, !expert, image) ;
308 Add2RecPointsList(h13, 13, expert, !image) ;
309 Add2RecPointsList(h14, 14, expert, !image) ;
310 Add2RecPointsList(h15, 15, expert, !image) ;
311 Add2RecPointsList(h16, 16, expert, !image) ;
312 Add2RecPointsList(h17, 17, expert, !image) ;
313 Add2RecPointsList(h18, 18, expert, !image) ;
314 Add2RecPointsList(h19, 19, expert, !image) ;
317 //____________________________________________________________________________
318 void AliTOFQADataMakerRec::InitESDs()
321 //create ESDs histograms in ESDs subdir
324 const Bool_t expert = kTRUE ;
325 const Bool_t image = kTRUE ;
327 TH1F * h0 = new TH1F("hTOFESDs", "Number of matched TOF tracks per event;Number of TOF matched ESD tracks;Counts", 201, -1., 200.) ;
328 TH1F * h1 = new TH1F("hTOFESDsTime", "Matched ESDs tracks: TOF Time spectrum; Calibrated TOF time [ns];Counts", 25000, 0., 610. ) ;
329 TH1F * h2 = new TH1F("hTOFESDsRawTime", "Matched ESDs tracks: TOF raw Time spectrum;Measured TOF time [ns];Counts", 25000, 0., 610.) ;
330 TH1F * h3 = new TH1F("hTOFESDsToT", "Matched ESDs tracks: TOF ToT spectrum; ESDs ToT [ns];Counts",1000, 0., 48.8) ;
331 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.) ;
332 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.) ;
333 TH1F * h6 = new TH1F("hTOFESDsTPCmatched", "TPC-TOF matched tracks momentum distribution (GeV/c); p (GeV/c);Counts", 50, 0.20, 5.00) ;
334 TH1F * h7 = new TH1F("hTOFESDsMatchingProb", "TPC-TOF track-matching probability per event;TPC-TOF track-matching probability (%) ;Counts",101, -1.0, 100) ;
335 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) ;
336 TH1F * h9 = new TH1F("hTOFESDsDiffTime", "ESDs t_{TOF}-t_{exp} spectrum in TOF (ps); t_{TOF}-t_{exp} [ps];Counts", 200, -2440., 2440.) ;
337 TH1F * h10 = new TH1F("hTOFHitsLength", "Matched ESDs tracks: Length Spectrum; Track length [cm];Counts", 800, 0., 800) ;
350 Add2ESDsList(h0, 0, !expert, image) ;
351 Add2ESDsList(h1, 1, !expert, image) ;
352 Add2ESDsList(h2, 2, expert, image) ;
353 Add2ESDsList(h3, 3, !expert, image) ;
354 Add2ESDsList(h4, 4, expert, image) ;
355 Add2ESDsList(h5, 5, expert, image) ;
356 Add2ESDsList(h6, 6, expert, image) ;
357 Add2ESDsList(h7, 7, expert, image) ;
358 Add2ESDsList(h8, 8, expert, image) ;
359 Add2ESDsList(h9, 9, !expert, image) ;
360 Add2ESDsList(h10, 10, !expert, image) ;
364 //____________________________________________________________________________
365 void AliTOFQADataMakerRec::MakeRaws(AliRawReader* rawReader)
368 // makes data from Raws
370 if (rawReader->GetType()==7) {
372 Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;//in ns
373 Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
376 Int_t ntof[5]; /* 0=tot, 1=IA, 2=OA, 3=IC, 4=OC*/
377 for (Int_t j=0;j<5;j++){ ntof[j]=0;}
379 Int_t equipmentID[5]; //(ddl, trm, chain,tdc,channel)
380 Int_t volumeID[5]; //(sector,plate,strip,padX,padZ)
381 Int_t volumeID2[5]; //(sector,plate,strip,padZ,padX) to use AliTOFGeometry::GetIndex()
382 Int_t out[5]; // out=(indexZ,indexPhi)
385 TClonesArray * clonesRawData;
386 AliTOFRawStream tofInput(rawReader);
388 //uncomment if needed to apply DeltaBC correction
389 //tofInput.ApplyBCCorrections(kTRUE);
391 for (Int_t iDDL = 0; iDDL < AliTOFGeometry::NDDL()*AliTOFGeometry::NSectors(); iDDL++){
394 tofInput.LoadRawDataBuffers(iDDL);
395 clonesRawData = (TClonesArray*)tofInput.GetRawData();
396 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
398 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
400 if (tofRawDatum->GetTOF()){
403 equipmentID[1]=tofRawDatum->GetTRM();
404 equipmentID[2]=tofRawDatum->GetTRMchain();
405 equipmentID[3]=tofRawDatum->GetTDC();
406 equipmentID[4]=tofRawDatum->GetTDCchannel();
408 if (CheckEquipID(equipmentID)){
409 tofInput.EquipmentId2VolumeId(iDDL,
410 tofRawDatum->GetTRM(),
411 tofRawDatum->GetTRMchain(),
412 tofRawDatum->GetTDC(),
413 tofRawDatum->GetTDCchannel(),
415 if (FilterSpare(equipmentID)) continue;
416 if (FilterLTMData(equipmentID)){ //counts LTM hits
417 if (tofRawDatum->GetTOT()) GetRawsData(15)->Fill(equipmentID[0]);
419 if (CheckVolumeID(volumeID)){
421 GetMapIndeces(volumeID,out);
422 volumeID2[0]=volumeID[0];
423 volumeID2[1]=volumeID[1];
424 volumeID2[2]=volumeID[2];
425 volumeID2[3]=volumeID[4];
426 volumeID2[4]=volumeID[3];
427 chIndex=AliTOFGeometry::GetIndex(volumeID2);
429 if (tofRawDatum->GetTOT()){
430 if (!(fCalibData->GetNoiseStatus(chIndex)==AliTOFChannelOnlineStatusArray::kTOFNoiseBad)
431 && (fCalibData->GetHWStatus(chIndex) == AliTOFChannelOnlineStatusArray::kTOFHWOk)) {//noise and enabled filter
432 ntof[0]++; //counter for tof hits
434 //fill global spectra for DQM plots
435 GetRawsData(5)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;//in ns
436 GetRawsData(10)->Fill( tofRawDatum->GetTOT()*tot2ns) ;//in ns
438 //fill side-related spectra for experts plots
439 if (volumeID2[0]>4 && volumeID2[0]<14){ //I side
440 if ((iDDL%4==0)|| (iDDL%4==1)){ //A side
442 GetRawsData(6)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;
443 GetRawsData(11)->Fill( tofRawDatum->GetTOT()*tot2ns) ;
445 if ((iDDL%4==2)|| (iDDL%4==3)){//C side
447 GetRawsData(8)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;
448 GetRawsData(13)->Fill( tofRawDatum->GetTOT()*tot2ns) ;
452 if (volumeID2[0]<5 || volumeID2[0]>13){ //O side
453 if ((iDDL%4==0)|| (iDDL%4==1)){ //A side
455 GetRawsData(7)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;
456 GetRawsData(12)->Fill( tofRawDatum->GetTOT()*tot2ns) ;
458 if ((iDDL%4==2)|| (iDDL%4==3)){//C side
460 GetRawsData(9)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;
461 GetRawsData(14)->Fill( tofRawDatum->GetTOT()*tot2ns) ;
466 GetRawsData(18)->Fill(chIndex);
468 Int_t trm= iDDL*10+(equipmentID[1]-3);
469 if (iDDL>=0 && iDDL<36) {
470 GetRawsData(16)->Fill(trm) ;
471 GetRawsData(20)->Fill(trm,tofRawDatum->GetTOF()*tdc2ns);
472 GetRawsData(22)->Fill(trm,tofRawDatum->GetTOT()*tot2ns);
474 if (iDDL>=36 && iDDL<72) {
475 GetRawsData(17)->Fill(trm) ;//in ns
476 GetRawsData(21)->Fill(trm,tofRawDatum->GetTOF()*tdc2ns);
477 GetRawsData(23)->Fill(trm,tofRawDatum->GetTOT()*tot2ns);
479 GetRawsData(24)->Fill(GetStripIndex(volumeID),tofRawDatum->GetTOF()*tdc2ns) ;
480 //GetRawsData(25)->Fill( out[0],out[1]) ;//raw map
484 if (!(fCalibData->GetNoiseStatus(chIndex) == AliTOFChannelOnlineStatusArray::kTOFNoiseBad)
485 && (fCalibData->GetHWStatus(chIndex) == AliTOFChannelOnlineStatusArray::kTOFHWOk))
486 GetRawsData(19)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;//in ns
488 }//end volumeID check
493 clonesRawData->Clear();
496 for (Int_t j=0;j<5;j++){
497 if(ntof[j]<=0) GetRawsData(j)->Fill(-1.) ;
498 else GetRawsData(j)->Fill(ntof[j]);
500 fProcessedRawEventN++;
503 AliDebug(1,Form("Event of type %d found. Skipping non-physics event for QA.\n", rawReader->GetType()));
505 EnableDqmShifterOpt(kTRUE);
508 //____________________________________________________________________________
509 void AliTOFQADataMakerRec::MakeRecPoints(TTree * clustersTree)
512 // Make data from Clusters
515 Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
516 Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
522 TBranch *branch=clustersTree->GetBranch("TOF");
524 AliError("can't get the branch with the TOF clusters !");
528 static TClonesArray dummy("AliTOFcluster",10000);
530 TClonesArray *clusters=&dummy;
531 branch->SetAddress(&clusters);
534 clustersTree->GetEvent(0);
536 Int_t nentries=clusters->GetEntriesFast();
538 GetRecPointsData(0)->Fill(-1.) ;
540 GetRecPointsData(0)->Fill(nentries) ;
543 TIter next(clusters) ;
545 while ( (c = dynamic_cast<AliTOFcluster *>(next())) ) {
547 volumeID[0] = c->GetDetInd(0);
548 volumeID[1] = c->GetDetInd(1);
549 volumeID[2] = c->GetDetInd(2);
550 volumeID[3] = c->GetDetInd(4); //X and Z indeces inverted in RecPoints
551 volumeID[4] = c->GetDetInd(3); //X and Z indeces inverted in RecPoints
553 for (Int_t i=0;i<5;i++){
554 volumeID2[i]=c->GetDetInd(i); //X and Z indeces inverted in RecPoints
557 GetMapIndeces(volumeID,out);
558 Int_t chIndex=AliTOFGeometry::GetIndex(volumeID2);
559 Int_t iDDL=AliTOFRawStream::Geant2DDL(volumeID2);
560 Int_t iTRM=AliTOFRawStream::Geant2TRM(volumeID2);
562 if (c->GetTDCRAW() && c->GetTDC() && c->GetToT()){
563 if (volumeID2[0]>4 && volumeID2[0]<14){ //I side
564 if ((iDDL%4==0)|| (iDDL%4==1)){ //A side
565 GetRecPointsData(1)->Fill( c->GetTDC()*tdc2ns) ;//in ns
566 GetRecPointsData(5)->Fill( c->GetTDCRAW()*tdc2ns) ;//in ns
567 GetRecPointsData(9)->Fill( c->GetToT()*tot2ns) ;//in ns
570 if ((iDDL%4==2)|| (iDDL%4==3)){//C side
571 GetRecPointsData(3)->Fill( c->GetTDC()*tdc2ns) ;//in ns
572 GetRecPointsData(7)->Fill( c->GetTDCRAW()*tdc2ns) ;//in ns
573 GetRecPointsData(11)->Fill( c->GetToT()*tot2ns) ;//in ns
577 if (volumeID2[0]<5 || volumeID2[0]>13){ //O side
578 if ((iDDL%4==0)|| (iDDL%4==1)){ //A side
579 GetRecPointsData(2)->Fill( c->GetTDC()*tdc2ns) ;//in ns
580 GetRecPointsData(6)->Fill( c->GetTDCRAW()*tdc2ns) ;//in ns
581 GetRecPointsData(10)->Fill( c->GetToT()*tot2ns) ;//in ns
583 if ((iDDL%4==2)|| (iDDL%4==3)){//C side
584 GetRecPointsData(4)->Fill( c->GetTDC()*tdc2ns) ;//in ns
585 GetRecPointsData(8)->Fill( c->GetTDCRAW()*tdc2ns) ;//in ns
586 GetRecPointsData(12)->Fill( c->GetToT()*tot2ns) ;//in ns
591 GetRecPointsData(14)->Fill(out[0],out[1]);
592 GetRecPointsData(15)->Fill(GetStripIndex(volumeID), c->GetTDC()*tdc2ns) ;
593 Int_t trm= iDDL*10+(iTRM-3);
594 if (iDDL>=0 && iDDL<36) {
595 GetRecPointsData(16)->Fill(trm,c->GetTDC()*tdc2ns);
596 GetRecPointsData(18)->Fill(trm,c->GetToT()*tot2ns);
598 if (iDDL>=36 && iDDL<72) {
599 GetRecPointsData(17)->Fill(trm,c->GetTDC()*tdc2ns);
600 GetRecPointsData(19)->Fill(trm,c->GetToT()*tot2ns);
602 GetRecPointsData(13)->Fill(chIndex) ;//in ns
606 EnableDqmShifterOpt(kFALSE);
609 //____________________________________________________________________________
610 void AliTOFQADataMakerRec::MakeESDs(AliESDEvent * esd)
613 // make QA data from ESDs
615 const Double_t speedOfLight = TMath::C()*1E2*1E-12; // cm/ps
616 const Double_t pionMass = 0.13957018; //GeV/c^2
618 Int_t ntrk = esd->GetNumberOfTracks() ;
627 AliESDtrack *track=esd->GetTrack(ntrk);
628 Double_t tofTime=track->GetTOFsignal();//in ps
629 Double_t tofTimeRaw=track->GetTOFsignalRaw();//in ps
630 Double_t tofToT=track->GetTOFsignalToT(); //in ps
632 UInt_t status=track->GetStatus();
633 if (track->IsOn(AliESDtrack::kTPCrefit)) {
635 Double_t y=track->Eta();
636 if ((status&AliESDtrack::kTOFout)!=0) {
641 if (TMath::Abs(y)<0.9) {
642 GetESDsData(6)->Fill(track->Pt());
645 GetESDsData(1)->Fill(tofTime*1E-3);
646 GetESDsData(2)->Fill(tofTimeRaw*1E-3);
647 GetESDsData(3)->Fill(tofToT*1E-3);
648 //check how many tracks where ESD PID is ok
649 if ((status&AliESDtrack::kESDpid)!=0) nesdpid++;
650 if (((status&AliESDtrack::kESDpid)&AliESDtrack::kTOFpid)!=0) ntofpid++;
652 Double_t length =track->GetIntegratedLength();
653 Double_t mom2=(track->Pt()*track->Pt())+(track->Pz()*track->Pz());
654 Double_t eTexp = TMath::Sqrt(1+(pionMass*pionMass/mom2))*length/speedOfLight; //in ps
655 GetESDsData(9)->Fill(tofTime-eTexp);
656 GetESDsData(10)->Fill(length);
657 } //end check on matched tracks
659 }//end check on TPCrefit
663 GetESDsData(0)->Fill(-1.) ;
665 GetESDsData(0)->Fill(nentries) ;
669 Float_t ratio = (Int_t)ntofpid/(Int_t)ntof*100; //identified by TOF over matched
670 GetESDsData(4)->Fill(ratio) ;
674 Float_t ratio = (Float_t)ntofpid/(Float_t)nesdpid *100; //identified by TOF over identified
675 GetESDsData(5)->Fill(ratio) ;
679 Float_t ratio = (Float_t)ntof/(Float_t)ntpc*100.; //matching probability
680 GetESDsData(7)->Fill(ratio) ;
684 Float_t ratio = (Float_t)ntof/(Float_t)ntofout*100; //matched over propagated to TOF outer radius
685 GetESDsData(8)->Fill(ratio) ;
687 EnableDqmShifterOpt(kFALSE);
690 //____________________________________________________________________________
691 void AliTOFQADataMakerRec::StartOfDetectorCycle()
694 //Detector specific actions at start of cycle
696 fCalibData = GetCalibData();
700 //____________________________________________________________________________
701 void AliTOFQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
703 //Detector specific actions at end of cycle
704 // do the QA checking
706 for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
707 if ( !AliQAv1::Instance()->IsEventSpecieSet(specie) )
710 AliInfo(Form("Processed %i physics raw events",fProcessedRawEventN));
712 if (fEnableDqmShifterOpt){
713 // Help make the raw qa histogram easier to interpret for the DQM shifter
714 // This is still to be optimized...
716 if (!GetRawsData(0) || !GetRawsData(5) || !GetRawsData(10)
717 || !GetRawsData(15) || !GetRawsData(16) || !GetRawsData(17)) {
718 printf("No histogram for DQM found - Possible memory corruption ???. Please check\n") ;
722 Double_t yTimeMin = GetRawsData(5)->GetMinimum();
723 Double_t yTimeMax = GetRawsData(5)->GetMaximum();
724 Double_t yTotMin = GetRawsData(10)->GetMinimum();
725 Double_t yTotMax = GetRawsData(10)->GetMaximum();
727 TLine* lineExpTimeMin = new TLine(200., yTimeMin, 200., yTimeMax);
728 lineExpTimeMin->SetLineColor(kGreen);
729 lineExpTimeMin->SetLineWidth(2);
731 TLine* lineExpTimeMax = new TLine(250., yTimeMin, 250., yTimeMax);
732 lineExpTimeMax->SetLineColor(kGreen);
733 lineExpTimeMax->SetLineWidth(2);
735 TLine* lineExpTotMin = new TLine( 5., yTotMin, 5., yTotMax);
736 lineExpTotMin->SetLineColor(kGreen);
737 lineExpTotMin->SetLineWidth(2);
739 TLine* lineExpTotMax = new TLine(20., yTotMin, 20., yTotMax);
740 lineExpTotMax->SetLineColor(kGreen);
741 lineExpTotMax->SetLineWidth(2);
743 for(Int_t j=0;j<5;j++){
744 //make up for hits count
745 GetRawsData(j)->SetMarkerColor(kRed);
746 GetRawsData(j)->SetMarkerStyle(7);
748 //make up for time spectra
749 ((TH1F*)GetRawsData(j+5))->GetListOfFunctions()->Add(lineExpTimeMin);
750 ((TH1F*)GetRawsData(j+5))->GetListOfFunctions()->Add(lineExpTimeMax);
751 GetRawsData(j+5)->SetLineColor(kBlue+1);
752 GetRawsData(j+5)->SetMarkerColor(kBlue+1);
754 //make up for tot spectra
755 ((TH1F*)GetRawsData(j+10))->GetListOfFunctions()->Add(lineExpTotMin);
756 ((TH1F*)GetRawsData(j+10))->GetListOfFunctions()->Add(lineExpTotMax);
757 GetRawsData(j+10)->SetLineColor(kBlue+1);
758 GetRawsData(j+10)->SetMarkerColor(kBlue+1);
761 //make up for equipment hits count
762 for(Int_t j=15;j<18;j++){
763 GetRawsData(j)->SetLineColor(kBlue+1);
764 GetRawsData(j)->SetLineWidth(1);
766 GetRawsData(j)->SetMarkerStyle(8);
767 GetRawsData(j)->SetMarkerSize(0.7);
768 GetRawsData(j)->SetMarkerColor(kBlue+2);
770 Int_t ySMmax=GetRawsData(j)->GetMaximum();
773 for (Int_t sm=0;sm<10;sm++){
774 lineSMid[sm] = new TLine( 40*sm+360*(j%16), 0, 40*sm+360*(j%16), ySMmax);
775 lineSMid[sm]->SetLineColor(kMagenta);
776 lineSMid[sm]->SetLineWidth(2);
777 GetRawsData(j)->GetListOfFunctions()->Add(lineSMid[sm]);
778 GetRawsData(j)->SetFillColor(kBlue+1);
783 //make up for orphans histo
784 GetRawsData(19)->SetLineColor(kBlue+1);
785 GetRawsData(19)->SetMarkerColor(kBlue+1);
786 ((TH1F*)GetRawsData(19))->GetListOfFunctions()->Add(lineExpTimeMin);
787 ((TH1F*)GetRawsData(19))->GetListOfFunctions()->Add(lineExpTimeMax);
790 AliQAChecker::Instance()->Run(AliQAv1::kTOF, task, list) ;
792 //____________________________________________________________________________
793 void AliTOFQADataMakerRec::GetMapIndeces(const Int_t* const in , Int_t* out)
796 //return appropriate indeces for the theta-phi map
799 Int_t npadX = AliTOFGeometry::NpadX();
800 Int_t npadZ = AliTOFGeometry::NpadZ();
801 Int_t nStripA = AliTOFGeometry::NStripA();
802 Int_t nStripB = AliTOFGeometry::NStripB();
803 Int_t nStripC = AliTOFGeometry::NStripC();
805 Int_t isector = in[0];
806 Int_t iplate = in[1];
807 Int_t istrip = in[2];
811 Int_t stripOffset = 0;
817 stripOffset = nStripC;
820 stripOffset = nStripC+nStripB;
823 stripOffset = nStripC+nStripB+nStripA;
826 stripOffset = nStripC+nStripB+nStripA+nStripB;
829 AliDebug(1,Form("Wrong plate number in TOF (%d) !",iplate));
832 Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
833 Int_t phiindex=npadX*isector+ipadX+1;
839 //---------------------------------------------------------------
840 Int_t AliTOFQADataMakerRec::GetStripIndex(const Int_t * const in)
842 /* return tof strip index between 0 and 91 */
844 Int_t nStripA = AliTOFGeometry::NStripA();
845 Int_t nStripB = AliTOFGeometry::NStripB();
846 Int_t nStripC = AliTOFGeometry::NStripC();
848 // Int_t isector = in[0];
849 Int_t iplate = in[1];
850 Int_t istrip = in[2];
851 //Int_t ipadX = in[3];
852 //Int_t ipadZ = in[4];
854 Int_t stripOffset = 0;
860 stripOffset = nStripC;
863 stripOffset = nStripC+nStripB;
866 stripOffset = nStripC+nStripB+nStripA;
869 stripOffset = nStripC+nStripB+nStripA+nStripB;
872 AliDebug(1,Form("Wrong plate number in TOF (%d) !",iplate));
877 if (stripOffset<0 || stripOffset>92) return -1;
879 return (stripOffset+istrip);
882 //---------------------------------------------------------------
883 Bool_t AliTOFQADataMakerRec::CheckVolumeID(const Int_t * const volumeID)
886 //Checks volume ID validity
889 for (Int_t j=0;j<5;j++){
891 AliDebug(1,Form("Invalid detector volume index for volumeID[%i]",j));
899 //---------------------------------------------------------------
900 Bool_t AliTOFQADataMakerRec::CheckEquipID(const Int_t * const equipmentID)
903 //Checks equipment ID validity
905 for (Int_t j=0;j<5;j++){
906 if (equipmentID[j]<0) {
907 AliDebug(1,Form("Invalid equipment volume index for equipmentID[%i]",j));
913 //---------------------------------------------------------------
914 Bool_t AliTOFQADataMakerRec::FilterLTMData(const Int_t * const equipmentID) const
916 /*It returns kTRUE if data come from LTM.
917 It thus filters trigger-related signals */
920 //if (!CheckEquipID(equipmentID)) return kFALSE;
921 ddl = equipmentID[0];
922 trm = equipmentID[1];
923 tdc = equipmentID[3];
925 if ((ddl%2==1) && (trm==3) && (tdc>11 && tdc<15))
931 //---------------------------------------------------------------
932 Bool_t AliTOFQADataMakerRec::FilterSpare(const Int_t * const equipmentID) const
934 /*It returns kTRUE if data come from spare
936 So far only check on TRM 3 crate left is implemented */
939 //if (!CheckEquipID(equipmentID)) return kFALSE;
940 ddl = equipmentID[0];
941 trm = equipmentID[1];
942 tdc = equipmentID[3];
944 if ((ddl%2==1) && (trm==3) && (tdc>2 && tdc<12))