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 // last modified by F. Bellini (fbellini@cern.ch) on 02/03/2010 //
23 ///////////////////////////////////////////////////////////////////////
25 /* Modified by fbellini on 30/03/2010
26 Changed raws time histos range to 610ns
27 Added FilterLTMData() and FilterSpare() methods
28 Added check on enabled channels for raw data
33 /* Modified by fbellini on 02/03/2010
34 Fixed raw data decoding methods (use AliTOFRawStream::LoadRawDataBuffer())
35 Added filter for noisy channels and read map from OCDB
36 Added GetCalibData() method
37 Added CheckVolumeID() and CheckEquipID() methods
41 #include <TClonesArray.h>
46 #include "AliCDBManager.h"
47 #include "AliCDBEntry.h"
48 #include "AliESDEvent.h"
49 #include "AliESDtrack.h"
50 #include "AliQAChecker.h"
51 #include "AliRawReader.h"
52 #include "AliTOFRawStream.h"
53 #include "AliTOFcluster.h"
54 #include "AliTOFQADataMakerRec.h"
55 #include "AliTOFrawData.h"
56 #include "AliTOFGeometry.h"
57 #include "AliTOFChannelOnlineStatusArray.h"
60 ClassImp(AliTOFQADataMakerRec)
62 //____________________________________________________________________________
63 AliTOFQADataMakerRec::AliTOFQADataMakerRec() :
64 AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kTOF), "TOF Quality Assurance Data Maker"),
66 fEnableNoiseFiltering(kFALSE)
73 //____________________________________________________________________________
74 AliTOFQADataMakerRec::AliTOFQADataMakerRec(const AliTOFQADataMakerRec& qadm) :
76 fCalibData(qadm.fCalibData),
77 fEnableNoiseFiltering(qadm.fEnableNoiseFiltering)
82 SetName((const char*)qadm.GetName()) ;
83 SetTitle((const char*)qadm.GetTitle());
86 //__________________________________________________________________
87 AliTOFQADataMakerRec& AliTOFQADataMakerRec::operator = (const AliTOFQADataMakerRec& qadm )
90 // assignment operator.
92 this->~AliTOFQADataMakerRec();
93 new(this) AliTOFQADataMakerRec(qadm);
97 //----------------------------------------------------------------------------
98 AliTOFChannelOnlineStatusArray* AliTOFQADataMakerRec::GetCalibData() const
101 // Retrive TOF calib objects from OCDB
103 AliCDBManager *man = AliCDBManager::Instance();
107 cdbe = man->Get("TOF/Calib/Status",fRun);
109 AliWarning("Load of calibration data from default storage failed!");
110 AliWarning("Calibration data will be loaded from local storage ($ALICE_ROOT)");
111 man->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
112 cdbe = man->Get("TOF/Calib/Status",fRun);
114 // Retrieval of data in directory TOF/Calib/Data:
116 AliTOFChannelOnlineStatusArray * array = 0;
117 if (cdbe) array = (AliTOFChannelOnlineStatusArray *)cdbe->GetObject();
118 if (!array) AliFatal("No calibration data from calibration database !");
126 //____________________________________________________________________________
127 void AliTOFQADataMakerRec::InitRaws()
130 // create Raws histograms in Raws subdir
132 const Bool_t expert = kTRUE ;
133 const Bool_t saveCorr = kTRUE ;
134 const Bool_t image = kTRUE ;
136 TH1F * h0 = new TH1F("hTOFRaws", "Number of TOF Raw Hits; TOF raw hits number;Counts ",1001, -1., 1000.) ;
137 TH1F * h1 = new TH1F("hTOFRawsIA", "Number of TOF Raw Hits - I/A side; TOF raw hits number ;Counts ",1001, -1., 1000.) ;
138 TH1F * h2 = new TH1F("hTOFRawsOA", "Number of TOF Raw Hits - O/A side; TOF raw hits number ;Counts ",1001, -1., 1000.) ;
139 TH1F * h3 = new TH1F("hTOFRawsIC", "Number of TOF Raw Hits - I/C side; TOF raw hits number ;Counts ",1001, -1., 1000.) ;
140 TH1F * h4 = new TH1F("hTOFRawsOC", "Number of TOF Raw Hits - O/C side; TOF raw hits number ;Counts ",1001, -1., 1000.) ;
141 TH1F * h5 = new TH1F("hTOFRawsTimeIA", "Raws Time Spectrum in TOF (ns) - I/A side;Measured raw TOF time [ns];Counts", 25000,0. ,610.) ;
142 TH1F * h6 = new TH1F("hTOFRawsTimeOA", "Raws Time Spectrum in TOF (ns) - O/A side;Measured raw TOF time [ns];Counts", 25000,0. ,610.) ;
143 TH1F * h7 = new TH1F("hTOFRawsTimeIC", "Raws Time Spectrum in TOF (ns) - I/C side;Measured raw TOF time [ns];Counts", 25000,0. ,610.) ;
144 TH1F * h8 = new TH1F("hTOFRawsTimeOC", "Raws Time Spectrum in TOF (ns) - O/C side;Measured raw TOF time [ns];Counts", 25000,0. ,610.) ;
145 TH1F * h9 = new TH1F("hTOFRawsToTIA", "Raws ToT Spectrum in TOF (ns) - I/A side;Measured raw ToT [ns];Counts", 1000, 0., 244.) ;
146 TH1F * h10 = new TH1F("hTOFRawsToTOA", "Raws ToT Spectrum in TOF (ns) - O/A side;Measured raw ToT [ns];Counts", 10000, 0., 244.) ;
147 TH1F * h11 = new TH1F("hTOFRawsToTIC", "Raws ToT Spectrum in TOF (ns) - I/C side;Measured raw ToT [ns];Counts", 10000, 0., 244.) ;
148 TH1F * h12 = new TH1F("hTOFRawsToTOC", "Raws ToT Spectrum in TOF (ns) - O/C side;Measured raw ToT [ns];Counts", 10000, 0., 244.) ;
150 TH1F * h13 = new TH1F("hTOFRawsTRMHitRate035", "TRM hit rate - crates 0 to 35 ;TRM index = DDL*10+TRM(0-9);Counts", 361, 0., 361.) ;
151 TH1F * h14 = new TH1F("hTOFRawsTRMHitRate3671","TRM hit rate - crates 36 to 71 ;TRM index = DDL*10+TRM(0-9);Counts", 361, 361., 722.) ;
152 TH1F * h15 = new TH1F("hTOFRawsLTMHitRate", "LTM hit rate; Crate; Counts", 72, 0., 72.);
154 TH1F * h16 = new TH1F("hTOFRawsVolumeErrorCounter","Invalid hit/volume association error counter per DDL; DDL; error counter ",73, -1., 72.) ;
155 TH1F * h17 = new TH1F("hTOFRawsEquipErrorCounter", "Invalid hit/equipment association error counter per DDL; DDL; error counter ",73, -1., 72.) ;
157 TH1F * h18 = new TH1F("hTOFOrphansTime", "Raws Time Spectrum in TOF (ns) for hits with no ToT measurement;Measured raw TOF time [ns];Counts", 25000,0. ,610.) ;
158 TH1F * h19 = new TH1F("hTOFRawsWords", "Number of TOF Raw words per event; TOF raw words number;Counts ",1001, -1., 1000.) ;
159 TH2F * h20 = new TH2F("hTOFRawTimeVsDDL", "TOF raw time spectrum vs DDL; TOF raw time [ns];DDL ",1000,0. ,2440., 72, 0., 72.) ;
160 TH2F * h21 = new TH2F("hTOFRawToTVsDDL", "TOF raw ToT spectrum vs DDL; TOF raw ToT [ns];DDL ", 100, 0., 244., 72, 0., 72.) ;
161 TH2F * h22 = 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) ;
163 // TH1F * h23 = new TH1F("hTOFRawsDecodingErr","Decoding errors counts per DDL; DDL;Counts",72,0,72);
164 // TH1F * h24 = new TH1F("hTOFRawsEventsPerDDL","Events number per DDL; DDL; Event Number",72,0,72);
191 Add2RawsList(h0, 0, !expert, image, !saveCorr) ;
192 Add2RawsList(h1, 1, expert, !image, !saveCorr) ;
193 Add2RawsList(h2, 2, expert, !image, !saveCorr) ;
194 Add2RawsList(h3, 3, expert, !image, !saveCorr) ;
195 Add2RawsList(h4, 4, expert, !image, !saveCorr) ;
196 Add2RawsList(h5, 5, !expert, image, !saveCorr) ;
197 Add2RawsList(h6, 6, !expert, image, !saveCorr) ;
198 Add2RawsList(h7, 7, !expert, image, !saveCorr) ;
199 Add2RawsList(h8, 8, !expert, image, !saveCorr) ;
200 Add2RawsList(h9, 9, !expert, image, !saveCorr) ;
201 Add2RawsList(h10, 10, !expert, image, !saveCorr) ;
202 Add2RawsList(h11, 11, !expert, image, !saveCorr) ;
203 Add2RawsList(h12, 12, !expert, image, !saveCorr) ;
204 Add2RawsList(h13, 13, expert, !image, !saveCorr) ;
205 Add2RawsList(h14, 14, expert, !image, !saveCorr) ;
206 Add2RawsList(h15, 15, expert, !image, !saveCorr) ;
207 Add2RawsList(h16, 16, expert, !image, !saveCorr) ;
208 Add2RawsList(h17, 17, expert, !image, !saveCorr) ;
209 Add2RawsList(h18, 18, expert, !image, !saveCorr) ;
210 Add2RawsList(h19, 19, expert, !image, !saveCorr) ;
211 Add2RawsList(h20, 20, expert, !image, !saveCorr) ;
212 Add2RawsList(h21, 21, expert, !image, !saveCorr) ;
213 Add2RawsList(h22, 22, !expert, image, !saveCorr) ;
218 //____________________________________________________________________________
219 void AliTOFQADataMakerRec::InitRecPoints()
222 // create RecPoints histograms in RecPoints subdir
225 const Bool_t expert = kTRUE ;
226 const Bool_t image = kTRUE ;
228 TH1F * h0 = new TH1F("hTOFRecPoints", "Number of TOF RecPoints per event;TOF recPoints number;Counts",1001, -1., 1000.) ;
229 TH1F * h1 = new TH1F("hTOFRecPointsTimeIA", "RecPoints Time Spectrum in TOF (ns) - I/A side;Calibrated TOF time [ns];Counts", 25000,0. ,610.) ;
230 TH1F * h2 = new TH1F("hTOFRecPointsTimeOA", "RecPoints Time Spectrum in TOF (ns)- O/A side;Calibrated TOF time [ns];Counts", 25000,0. ,610.) ;
231 TH1F * h3 = new TH1F("hTOFRecPointsTimeIC", "RecPoints Time Spectrum in TOF (ns)- I/C side;Calibrated TOF time [ns];Counts", 25000,0. ,610.) ;
232 TH1F * h4 = new TH1F("hTOFRecPointsTimeOC", "RecPoints Time Spectrum in TOF (ns)- O/C side;Calibrated TOF time [ns];Counts", 25000,0. ,610.) ;
234 TH1F * h5 = new TH1F("hTOFRecPointsRawTimeIA", "RecPoints raw Time Spectrum in TOF (ns)-I/A side ;Measured TOF time [ns];Counts", 25000,0. ,610.) ;
235 TH1F * h6 = new TH1F("hTOFRecPointsRawTimeOA", "RecPoints raw Time Spectrum in TOF (ns)-O/A side;Measured TOF time [ns];Counts", 25000,0. ,610.) ;
236 TH1F * h7 = new TH1F("hTOFRecPointsRawTimeIC", "RecPoints raw Time Spectrum in TOF (ns)-I/C side;Measured TOF time [ns];Counts", 25000,0. ,610.) ;
237 TH1F * h8 = new TH1F("hTOFRecPointsRawTimeOC", "RecPoints raw Time Spectrum in TOF (ns)-O/C side;Measured TOF time [ns];Counts", 25000,0. ,610.) ;
239 TH1F * h9 = new TH1F("hTOFRecPointsToTIA", "RecPoints ToT Spectrum in TOF (ns)-I/A side;Measured TOT [ns];Counts",10000, 0., 244. ) ;
240 TH1F * h10 = new TH1F("hTOFRecPointsToTOA", "RecPoints ToT Spectrum in TOF (ns)-O/A side;Measured TOT [ns];Counts",10000, 0., 244. ) ;
241 TH1F * h11 = new TH1F("hTOFRecPointsToTIC", "RecPoints ToT Spectrum in TOF (ns)-I/C side;Measured TOT [ns];Counts",10000, 0., 244. ) ;
242 TH1F * h12 = new TH1F("hTOFRecPointsToTOC", "RecPoints ToT Spectrum in TOF (ns)-O/C side;Measured TOT [ns];Counts",10000, 0., 244. ) ;
244 TH2F * h13 = 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) ;
261 Add2RecPointsList(h0, 0, !expert, image) ;
262 Add2RecPointsList(h1, 1, !expert, image) ;
263 Add2RecPointsList(h2, 2, !expert, image) ;
264 Add2RecPointsList(h3, 3, !expert, image) ;
265 Add2RecPointsList(h4, 4, !expert, image) ;
266 Add2RecPointsList(h5, 5, expert, !image) ;
267 Add2RecPointsList(h6, 6, expert, !image) ;
268 Add2RecPointsList(h7, 7, expert, !image) ;
269 Add2RecPointsList(h8, 8, expert, !image) ;
270 Add2RecPointsList(h9, 9, !expert, image) ;
271 Add2RecPointsList(h10, 10, !expert, image) ;
272 Add2RecPointsList(h11, 11, !expert, image) ;
273 Add2RecPointsList(h12, 12, !expert, image) ;
274 Add2RecPointsList(h13, 13, expert, !image) ;
278 //____________________________________________________________________________
279 void AliTOFQADataMakerRec::InitESDs()
282 //create ESDs histograms in ESDs subdir
285 Bool_t expert = kFALSE;
287 TH1F * h0 = new TH1F("hTOFESDs", "Number of matched TOF tracks over ESDs", 250, -1., 4.) ;
289 Add2ESDsList(h0, 0, expert) ;
291 TH1F * h1 = new TH1F("hTOFESDsTime", "Time Spectrum in TOF (ns)", 1000, 0., 244) ;
293 Add2ESDsList(h1, 1, expert) ;
295 TH1F * h2 = new TH1F("hTOFESDsRawTime", "raw Time Spectrum in TOF (ns)", 1000, 0., 244) ;
297 Add2ESDsList(h2, 2, expert) ;
299 TH1F * h3 = new TH1F("hTOFESDsToT", "ToT Spectrum in TOF (ns)", 1000, 0., 244) ;
301 Add2ESDsList(h3, 3, expert) ;
303 TH1F * h4 = new TH1F("hTOFESDsPID", "Fraction of matched TOF tracks with good PID glag", 100, 0., 1.) ;
305 Add2ESDsList(h4, 4, expert) ;
308 //____________________________________________________________________________
309 void AliTOFQADataMakerRec::MakeRaws(AliRawReader* rawReader)
312 // makes data from Raws
315 Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
316 Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
319 Int_t ntof[5]; /* 0=tot, 1=IA, 2=OA, 3=IC, 4=OC*/
320 for (Int_t j=0;j<5;j++){ ntof[j]=0;}
322 Int_t equipmentID[5]; //(ddl, trm, chain,tdc,channel)
323 Int_t volumeID[5]; //(sector,plate,strip,padX,padZ)
324 Int_t volumeID2[5]; //(sector,plate,strip,padZ,padX) to use AliTOFGeometry::GetIndex()
325 Int_t out[5]; // out=(indexZ,indexPhi)
328 TClonesArray * clonesRawData;
329 AliTOFRawStream tofInput(rawReader);
331 //uncomment if needed to apply DeltaBC correction
332 //tofInput.ApplyBCCorrections(kTRUE);
334 for (Int_t iDDL = 0; iDDL < AliTOFGeometry::NDDL()*AliTOFGeometry::NSectors(); iDDL++){
336 tofInput.LoadRawDataBuffers(iDDL);
337 clonesRawData = (TClonesArray*)tofInput.GetRawData();
338 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
340 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
342 if (tofRawDatum->GetTOF()){
345 equipmentID[1]=tofRawDatum->GetTRM();
346 equipmentID[2]=tofRawDatum->GetTRMchain();
347 equipmentID[3]=tofRawDatum->GetTDC();
348 equipmentID[4]=tofRawDatum->GetTDCchannel();
350 if (!CheckEquipID(equipmentID)) GetRawsData(17)->Fill(equipmentID[0]);
352 tofInput.EquipmentId2VolumeId(iDDL,
353 tofRawDatum->GetTRM(),
354 tofRawDatum->GetTRMchain(),
355 tofRawDatum->GetTDC(),
356 tofRawDatum->GetTDCchannel(),
358 if (FilterSpare(equipmentID)) continue;
359 if (FilterLTMData(equipmentID)) GetRawsData(15)->Fill(equipmentID[0]);
361 if (!CheckVolumeID(volumeID)) GetRawsData(16)->Fill(equipmentID[0]);
364 GetMapIndeces(volumeID,out);
366 volumeID2[0]=volumeID[0];
367 volumeID2[1]=volumeID[1];
368 volumeID2[2]=volumeID[2];
369 volumeID2[3]=volumeID[4];
370 volumeID2[4]=volumeID[3];
371 chIndex=AliTOFGeometry::GetIndex(volumeID2);
373 if (tofRawDatum->GetTOT()){
374 if (!(fCalibData->GetNoiseStatus(chIndex)==AliTOFChannelOnlineStatusArray::kTOFNoiseBad)
375 && (fCalibData->GetHWStatus(chIndex) == AliTOFChannelOnlineStatusArray::kTOFHWOk)) {//noise and enabled filter
376 ntof[0]++; //counter for tof hits
377 if (volumeID2[0]>4 && volumeID2[0]<14){ //I side
378 if ((iDDL%4==0)|| (iDDL%4==1)){ //A side
380 GetRawsData(5)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;//in ns
381 GetRawsData(9)->Fill( tofRawDatum->GetTOT()*tot2ns) ;//in ns
383 if ((iDDL%4==2)|| (iDDL%4==3)){//C side
385 GetRawsData(7)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;//in ns
386 GetRawsData(11)->Fill( tofRawDatum->GetTOT()*tot2ns) ;//in ns
390 if (volumeID2[0]<5 || volumeID2[0]>13){ //O side
391 if ((iDDL%4==0)|| (iDDL%4==1)){ //A side
393 GetRawsData(6)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;//in ns
394 GetRawsData(10)->Fill( tofRawDatum->GetTOT()*tot2ns) ;//in ns
396 if ((iDDL%4==2)|| (iDDL%4==3)){//C side
398 GetRawsData(8)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;//in ns
399 GetRawsData(12)->Fill( tofRawDatum->GetTOT()*tot2ns) ;//in ns
405 Int_t trm= iDDL*10+(equipmentID[1]-3);
406 if (iDDL>=0 && iDDL<36) GetRawsData(13)->Fill(trm) ;//in ns
407 if (iDDL>=36 && iDDL<72) GetRawsData(14)->Fill(trm) ;//in ns
408 GetRawsData(22)->Fill( out[0],out[1]) ;//raw map
410 GetRawsData(20)->Fill(tofRawDatum->GetTOF()*tdc2ns,iDDL);
411 GetRawsData(21)->Fill(tofRawDatum->GetTOT()*tot2ns,iDDL);
416 if (!(fCalibData->GetNoiseStatus(chIndex) == AliTOFChannelOnlineStatusArray::kTOFNoiseBad)
417 && (fCalibData->GetHWStatus(chIndex) == AliTOFChannelOnlineStatusArray::kTOFHWOk))
418 GetRawsData(18)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;//in ns
420 }//end volumeID check
426 clonesRawData->Clear();
429 for (Int_t j=0;j<5;j++){
430 if(ntof[j]<=0) GetRawsData(j)->Fill(-1.) ;
431 else GetRawsData(j)->Fill(ntof[j]);
433 if (nentries<=0) GetRawsData(19)->Fill(-1.) ;
434 else GetRawsData(19)->Fill(nentries) ;
438 //____________________________________________________________________________
439 void AliTOFQADataMakerRec::MakeRecPoints(TTree * clustersTree)
442 // Make data from Clusters
445 Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
446 Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
452 TBranch *branch=clustersTree->GetBranch("TOF");
454 AliError("can't get the branch with the TOF clusters !");
458 static TClonesArray dummy("AliTOFcluster",10000);
460 TClonesArray *clusters=&dummy;
461 branch->SetAddress(&clusters);
464 clustersTree->GetEvent(0);
466 Int_t nentries=clusters->GetEntriesFast();
468 GetRecPointsData(0)->Fill(-1.) ;
470 GetRecPointsData(0)->Fill(nentries) ;
473 TIter next(clusters) ;
475 while ( (c = dynamic_cast<AliTOFcluster *>(next())) ) {
477 volumeID[0] = c->GetDetInd(0);
478 volumeID[1] = c->GetDetInd(1);
479 volumeID[2] = c->GetDetInd(2);
480 volumeID[3] = c->GetDetInd(4); //X and Z indeces inverted in RecPoints
481 volumeID[4] = c->GetDetInd(3); //X and Z indeces inverted in RecPoints
483 for (Int_t i=0;i<5;i++){
484 volumeID2[i]=c->GetDetInd(i); //X and Z indeces inverted in RecPoints
487 GetMapIndeces(volumeID,out);
488 //Int_t chIndex=AliTOFGeometry::GetIndex(volumeID2);
489 Int_t iDDL=AliTOFRawStream::Geant2DDL(volumeID2);
491 if (c->GetTDCRAW() && c->GetTDC() && c->GetToT()){
492 if (volumeID2[0]>4 && volumeID2[0]<14){ //I side
493 if ((iDDL%4==0)|| (iDDL%4==1)){ //A side
494 GetRecPointsData(1)->Fill( c->GetTDC()*tdc2ns) ;//in ns
495 GetRecPointsData(5)->Fill( c->GetTDCRAW()*tdc2ns) ;//in ns
496 GetRecPointsData(9)->Fill( c->GetToT()*tot2ns) ;//in ns
499 if ((iDDL%4==2)|| (iDDL%4==3)){//C side
500 GetRecPointsData(3)->Fill( c->GetTDC()*tdc2ns) ;//in ns
501 GetRecPointsData(7)->Fill( c->GetTDCRAW()*tdc2ns) ;//in ns
502 GetRecPointsData(11)->Fill( c->GetToT()*tot2ns) ;//in ns
506 if (volumeID2[0]<5 || volumeID2[0]>13){ //O side
507 if ((iDDL%4==0)|| (iDDL%4==1)){ //A side
508 GetRecPointsData(2)->Fill( c->GetTDC()*tdc2ns) ;//in ns
509 GetRecPointsData(6)->Fill( c->GetTDCRAW()*tdc2ns) ;//in ns
510 GetRecPointsData(10)->Fill( c->GetToT()*tot2ns) ;//in ns
512 if ((iDDL%4==2)|| (iDDL%4==3)){//C side
513 GetRecPointsData(4)->Fill( c->GetTDC()*tdc2ns) ;//in ns
514 GetRecPointsData(8)->Fill( c->GetTDCRAW()*tdc2ns) ;//in ns
515 GetRecPointsData(12)->Fill( c->GetToT()*tot2ns) ;//in ns
520 GetRecPointsData(13)->Fill(out[0],out[1]);
529 //____________________________________________________________________________
530 void AliTOFQADataMakerRec::MakeESDs(AliESDEvent * const esd)
533 // make QA data from ESDs
535 Int_t ntrk = esd->GetNumberOfTracks() ;
539 AliESDtrack *track=esd->GetTrack(ntrk);
540 Double_t tofTime=track->GetTOFsignal()*1E-3;//in ns
541 Double_t tofTimeRaw=track->GetTOFsignalRaw()*1E-3;//in ns
542 Double_t tofToT=track->GetTOFsignalToT(); //in ns
543 if(!(tofTime>0))continue;
545 GetESDsData(1)->Fill(tofTime);
546 GetESDsData(2)->Fill(tofTimeRaw);
547 GetESDsData(3)->Fill(tofToT);
548 //check how many tracks where ESD PID is ok
549 UInt_t status=track->GetStatus();
550 if (((status&AliESDtrack::kESDpid)==0) ||
551 ((status&AliESDtrack::kTOFpid)==0)) continue;
557 GetESDsData(0)->Fill(-1.) ;
559 GetESDsData(0)->Fill(TMath::Log10(nentries)) ;
562 if(ntof>0)GetESDsData(4)->Fill(ntofpid/ntof) ;
565 //____________________________________________________________________________
566 void AliTOFQADataMakerRec::StartOfDetectorCycle()
569 //Detector specific actions at start of cycle
572 fCalibData = GetCalibData();
576 //____________________________________________________________________________
577 void AliTOFQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
579 //Detector specific actions at end of cycle
580 // do the QA checking
582 AliQAChecker::Instance()->Run(AliQAv1::kTOF, task, list) ;
584 //____________________________________________________________________________
585 void AliTOFQADataMakerRec::GetMapIndeces(const Int_t* const in , Int_t* out)
588 //return appropriate indeces for the theta-phi map
591 Int_t npadX = AliTOFGeometry::NpadX();
592 Int_t npadZ = AliTOFGeometry::NpadZ();
593 Int_t nStripA = AliTOFGeometry::NStripA();
594 Int_t nStripB = AliTOFGeometry::NStripB();
595 Int_t nStripC = AliTOFGeometry::NStripC();
597 Int_t isector = in[0];
598 Int_t iplate = in[1];
599 Int_t istrip = in[2];
603 Int_t stripOffset = 0;
609 stripOffset = nStripC;
612 stripOffset = nStripC+nStripB;
615 stripOffset = nStripC+nStripB+nStripA;
618 stripOffset = nStripC+nStripB+nStripA+nStripB;
621 AliDebug(1,Form("Wrong plate number in TOF (%d) !",iplate));
624 Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
625 Int_t phiindex=npadX*isector+ipadX+1;
630 //---------------------------------------------------------------
631 Bool_t AliTOFQADataMakerRec::CheckVolumeID(const Int_t * const volumeID)
634 //Checks volume ID validity
637 for (Int_t j=0;j<5;j++){
639 AliDebug(1,Form("Invalid detector volume index for volumeID[%i]",j));
647 //---------------------------------------------------------------
648 Bool_t AliTOFQADataMakerRec::CheckEquipID(const Int_t * const equipmentID)
651 //Checks equipment ID validity
653 for (Int_t j=0;j<5;j++){
654 if (equipmentID[j]<0) {
655 AliDebug(1,Form("Invalid equipment volume index for equipmentID[%i]",j));
661 //---------------------------------------------------------------
662 Bool_t AliTOFQADataMakerRec::FilterLTMData(const Int_t * const equipmentID) const
664 /*It returns kTRUE if data come from LTM.
665 It thus filter trigger-related signals */
668 //if (!CheckEquipID(equipmentID)) return kFALSE;
669 ddl = equipmentID[0];
670 trm = equipmentID[1];
671 tdc = equipmentID[3];
673 if ((ddl%2==1) && (trm==3) && (tdc>11 && tdc<15))
679 //---------------------------------------------------------------
680 Bool_t AliTOFQADataMakerRec::FilterSpare(const Int_t * const equipmentID) const
682 /*It returns kTRUE if data come from spare
684 So far only check on TRM 3 crate left is implemented */
687 //if (!CheckEquipID(equipmentID)) return kFALSE;
688 ddl = equipmentID[0];
689 trm = equipmentID[1];
690 tdc = equipmentID[3];
692 if ((ddl%2==1) && (trm==3) && (tdc>2 && tdc<12))