]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFQADataMakerRec.cxx
Changed scope for AliTOFRawStream variable. It solves the problem with some TClonesAr...
[u/mrichter/AliRoot.git] / TOF / AliTOFQADataMakerRec.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16 ///////////////////////////////////////////////////////////////////////
17 //                                                                   //
18 //  Produces the data needed to calculate the TOF quality assurance. //
19 //  QA objects are 1 & 2 Dimensional histograms.                     //
20 //  author: S.Arcelli                                                //
21 //                                                                   //
22 ///////////////////////////////////////////////////////////////////////
23
24 /*
25     Modified by fbellini on 13/09/2010
26   - Set TLines as private members
27   - Set image flag for expert histos
28
29   Modified by fbellini on 14/06/2010
30   - Updated plots
31   - use LoadRawDataBuffersV2()
32
33   Modified by fbellini on 10/05/2010
34   - Fixed EndOfDetectorCycle() memory corruption bug
35
36   Modified by fbellini on 22/04/2010
37    - Added filter for physics events
38
39    Modified by fbellini on 16/04/2010
40    - Added EnableDqmShifterOpt() 
41    - Modified EndOfDetectorCycle() with options for DQM         
42    - Updated ESDs QA
43
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
49
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  
55    - Updated Raw QA
56 */
57
58 #include <TClonesArray.h>
59 #include <TH1F.h> 
60 #include <TH2F.h> 
61 #include <TLine.h>
62
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"
75
76
77 ClassImp(AliTOFQADataMakerRec)
78            
79 //____________________________________________________________________________ 
80   AliTOFQADataMakerRec::AliTOFQADataMakerRec() : 
81   AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kTOF), "TOF Quality Assurance Data Maker"),
82   fCalibData(0x0),
83   fEnableNoiseFiltering(kFALSE),
84     fEnableDqmShifterOpt(kFALSE),
85     fProcessedRawEventN(0),
86     fLineExpTimeMin(0x0),
87     fLineExpTimeMax(0x0),
88     fLineExpTotMin(0x0),
89     fLineExpTotMax(0x0),
90   fTOFRawStream(AliTOFRawStream())
91 {
92   //
93   // ctor
94   //
95    
96   for (Int_t sm=0;sm<10;sm++){
97       fLineSMid035[sm]=0x0;
98       fLineSMid3671[sm]=0x0;
99   }
100     
101 }
102
103 //____________________________________________________________________________ 
104 AliTOFQADataMakerRec::AliTOFQADataMakerRec(const AliTOFQADataMakerRec& qadm) :
105   AliQADataMakerRec(),
106   fCalibData(qadm.fCalibData),
107   fEnableNoiseFiltering(qadm.fEnableNoiseFiltering),
108   fEnableDqmShifterOpt(qadm.fEnableDqmShifterOpt),
109   fProcessedRawEventN(qadm.fProcessedRawEventN),
110   fLineExpTimeMin(qadm.fLineExpTimeMin),
111   fLineExpTimeMax(qadm.fLineExpTimeMax),
112   fLineExpTotMin(qadm.fLineExpTotMin),
113   fLineExpTotMax(qadm.fLineExpTotMax),
114   fTOFRawStream(qadm.fTOFRawStream)
115 {
116   //
117   //copy ctor 
118   //
119   SetName((const char*)qadm.GetName()) ; 
120   SetTitle((const char*)qadm.GetTitle()); 
121    
122   for (Int_t sm=0;sm<10;sm++){
123       fLineSMid035[sm]=qadm.fLineSMid035[sm];
124       fLineSMid3671[sm]=qadm.fLineSMid3671[sm];
125   }
126 }
127
128 //__________________________________________________________________
129 AliTOFQADataMakerRec& AliTOFQADataMakerRec::operator = (const AliTOFQADataMakerRec& qadm )
130 {
131   //
132   // assignment operator.
133   //
134   this->~AliTOFQADataMakerRec();
135   new(this) AliTOFQADataMakerRec(qadm);
136   return *this;
137 }
138  
139 //----------------------------------------------------------------------------
140 AliTOFQADataMakerRec::~AliTOFQADataMakerRec()
141 {
142
143   fTOFRawStream.Clear();
144
145 }
146 //----------------------------------------------------------------------------
147 AliTOFChannelOnlineStatusArray* AliTOFQADataMakerRec::GetCalibData() const
148 {
149   //
150   // Retrive TOF calib objects from OCDB
151   //
152     AliCDBManager *man = AliCDBManager::Instance();
153
154     AliCDBEntry *cdbe=0;
155
156     cdbe = man->Get("TOF/Calib/Status",fRun);
157     if(!cdbe){
158         AliWarning("Load of calibration data from default storage failed!");
159         AliWarning("Calibration data will be loaded from local storage ($ALICE_ROOT)");
160         man->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
161         cdbe = man->Get("TOF/Calib/Status",fRun);
162     }
163     // Retrieval of data in directory TOF/Calib/Data:
164     
165     AliTOFChannelOnlineStatusArray * array = 0;
166     if (cdbe) array = (AliTOFChannelOnlineStatusArray *)cdbe->GetObject();
167     if (!array)  AliFatal("No calibration data from calibration database !");
168     
169     return array;
170     
171 }
172
173
174
175 //____________________________________________________________________________ 
176 void AliTOFQADataMakerRec::InitRaws()
177 {
178   //
179   // create Raws histograms in Raws subdir
180   //
181   const Bool_t expert   = kTRUE ; 
182   const Bool_t saveCorr = kTRUE ; 
183   const Bool_t image    = kTRUE ; 
184
185   TH1I * h0 =  new TH1I("hTOFRaws",      "TOF raw hit multiplicity; TOF raw hits number;Counts ",200, 0, 200) ;
186
187   TH1I * h1 =  new TH1I("hTOFRawsIA",    "TOF raw hit multiplicity - I/A side; TOF raw hits number ;Counts ",200, 0, 200) ;
188   TH1I * h2 =  new TH1I("hTOFRawsOA",    "TOF raw hit multiplicity - O/A side; TOF raw hits number ;Counts ",200, 0, 200) ;
189   TH1I * h3 =  new TH1I("hTOFRawsIC",    "TOF raw hit multiplicity - I/C side; TOF raw hits number ;Counts ",200, 0, 200) ;
190   TH1I * h4 =  new TH1I("hTOFRawsOC",    "TOF raw hit multiplicity - O/C side; TOF raw hits number ;Counts ",200, 0, 200) ;
191
192   TH1F * h5  = new TH1F("hTOFRawsTime", "TOF Raws - Hit time (ns);Measured Hit time [ns];Counts", 25000,0. ,610.) ; 
193
194   TH1F * h6  = new TH1F("hTOFRawsTimeIA", "TOF Raws - Hit time (ns) - I/A side;Measured Hit time [ns];Counts", 25000,0. ,610.) ; 
195   TH1F * h7  = new TH1F("hTOFRawsTimeOA", "TOF Raws - Hit time (ns) - O/A side;Measured Hit time [ns];Counts", 25000,0. ,610.) ; 
196   TH1F * h8  = new TH1F("hTOFRawsTimeIC", "TOF Raws - Hit time (ns) - I/C side;Measured Hit time [ns];Counts", 25000,0. ,610.) ; 
197   TH1F * h9  = new TH1F("hTOFRawsTimeOC", "TOF Raws - Hit time (ns) - O/C side;Measured Hit time [ns];Counts", 25000,0. ,610.) ; 
198
199   TH1F * h10  = new TH1F("hTOFRawsToT", "TOF Raws - Hit ToT (ns);Measured Hit ToT (ns);Counts", 1000, 0., 48.8) ; 
200
201   TH1F * h11  = new TH1F("hTOFRawsToTIA", "TOF Raws - Hit ToT (ns) - I/A side;Measured Hit ToT (ns);Counts", 1000, 0., 48.8) ; 
202   TH1F * h12  = new TH1F("hTOFRawsToTOA", "TOF Raws - Hit ToT (ns) - O/A side;Measured Hit ToT (ns);Counts", 1000, 0., 48.8) ; 
203   TH1F * h13  = new TH1F("hTOFRawsToTIC", "TOF Raws - Hit ToT (ns) - I/C side;Measured Hit ToT (ns);Counts", 1000, 0., 48.8) ; 
204   TH1F * h14  = new TH1F("hTOFRawsToTOC", "TOF Raws - Hit ToT (ns) - O/C side;Measured Hit ToT (ns);Counts", 1000, 0., 48.8) ; 
205
206   TH1F * h15 = new TH1F("hTOFRawsLTMHits", "LTMs OR signals; Crate; Counts",  72, 0., 72.);
207   TH1F * h16  = new TH1F("hTOFRawsTRMHits035", "TRM hits  - crates 0 to 35 ;TRM index = SMid(crate*10)+TRM(0-9);Hits",  361, 0., 361.) ;
208   TH1F * h17  = new TH1F("hTOFRawsTRMHits3671","TRM hits  - crates 36 to 71 ;TRM index = SMid(crate*10)+TRM(0-9);Hits", 361, 360., 721.) ;
209   
210   TH1F * h18 = new TH1F("hTOFRawChannelHits","TOF channel hits count; Channel ID; Hits", 158000, 0., 158000);
211   
212   TH1F * h19  = new TH1F("hTOFOrphansTime", "TOF Raws - Orphans time (ns);Measured Hit time [ns];Counts", 25000, 0. ,610.) ; 
213   TH2F * h20 = new TH2F("hTOFRawTimeVsTRM035", "TOF raws - Hit time vs TRM - crates 0 to 35; TRM index = DDL*10+TRM(0-9);TOF raw time [ns]", 361, 0., 361., 250, 0., 610.0) ;
214   TH2F * h21 = new TH2F("hTOFRawTimeVsTRM3671", "TOF raws - Hit time vs TRM - crates 36 to 72; TRM index = DDL**10+TRM(0-9);TOF raw time [ns]", 361, 360., 721., 250, 0., 610.0) ;
215   TH2F * h22 = new TH2F("hTOFRawToTVsTRM035",  "TOF raws - Hit ToT vs TRM - crates 0 to 35; TRM index = DDL*10+TRM(0-9);TOF raw ToT [ns] ", 361, 0., 361., 100, 0., 48.8) ;
216   TH2F * h23 = new TH2F("hTOFRawToTVsTRM3671",  "TOF raws - Hit ToT vs TRM - crates 36 to 72; TRM index = DDL*10+TRM(0-9);TOF raw ToT [ns] ", 361, 360., 721., 100, 0., 48.8) ;
217   TH2F * h24 = new TH2F("hTOFTimeVsStrip","TOF Raws - Hit time vs. strip (theta); Strip index;Raws TOF time (ns) ", 91,0.,91, 250, 0., 610.) ; 
218   
219 // TH2F * h25 = new TH2F("hTOFRawsClusMap","Raws vs TOF eta-phi;eta (2*strip+padz);phi (48*sector+padx)",183, -0.5, 182.5,865,-0.5,864.5) ; 
220
221   h0->Sumw2() ;
222   h1->Sumw2() ;
223   h2->Sumw2() ;
224   h3->Sumw2() ;
225   h4->Sumw2() ;
226   h5->Sumw2() ;
227   h6->Sumw2() ;
228   h7->Sumw2() ;
229   h8->Sumw2() ;
230   h9->Sumw2() ;
231   h10->Sumw2() ;
232   h11->Sumw2() ;
233   h12->Sumw2() ;
234   h13->Sumw2() ;
235   h14->Sumw2() ;
236   h15->Sumw2() ;
237   h16->Sumw2() ;
238   h17->Sumw2() ;
239   h18->Sumw2() ;
240   h19->Sumw2() ;
241   h20->Sumw2() ;
242   h21->Sumw2() ;
243   h22->Sumw2() ;
244   h23->Sumw2() ;
245   h24->Sumw2() ;
246
247   Add2RawsList(h0,   0, !expert,  image, !saveCorr) ;
248   Add2RawsList(h1,   1,  expert,  image, !saveCorr) ;
249   Add2RawsList(h2,   2,  expert,  image, !saveCorr) ;
250   Add2RawsList(h3,   3,  expert,  image, !saveCorr) ;
251   Add2RawsList(h4,   4,  expert,  image, !saveCorr) ;
252   Add2RawsList(h5,   5, !expert,  image, !saveCorr) ;
253   Add2RawsList(h6,   6,  expert,  image, !saveCorr) ;
254   Add2RawsList(h7,   7,  expert,  image, !saveCorr) ;
255   Add2RawsList(h8,   8,  expert,  image, !saveCorr) ;
256   Add2RawsList(h9,   9,  expert,  image, !saveCorr) ;
257   Add2RawsList(h10, 10, !expert,  image, !saveCorr) ;
258   Add2RawsList(h11, 11,  expert, !image, !saveCorr) ;
259   Add2RawsList(h12, 12,  expert, !image, !saveCorr) ;
260   Add2RawsList(h13, 13,  expert, !image, !saveCorr) ;
261   Add2RawsList(h14, 14,  expert, !image, !saveCorr) ;
262   Add2RawsList(h15, 15, !expert,  image, !saveCorr) ;
263   Add2RawsList(h16, 16, !expert,  image, !saveCorr) ;
264   Add2RawsList(h17, 17, !expert,  image, !saveCorr) ;
265   Add2RawsList(h18, 18,  expert, !image, !saveCorr) ;
266   Add2RawsList(h19, 19,  expert, !image, !saveCorr) ;
267   Add2RawsList(h20, 20,  expert, !image, !saveCorr) ;
268   Add2RawsList(h21, 21,  expert, !image, !saveCorr) ;
269   Add2RawsList(h22, 22,  expert, !image, !saveCorr) ;
270   Add2RawsList(h23, 23,  expert, !image, !saveCorr) ;
271   Add2RawsList(h24, 24,  expert, !image, !saveCorr) ;
272
273   //add lines for DQM shifter
274   fLineExpTimeMin = new TLine(200., 0., 200., 0.);
275   fLineExpTimeMin->SetLineColor(kGreen);
276   fLineExpTimeMin->SetLineWidth(2);
277   
278   fLineExpTimeMax = new TLine(250., 0., 250., 0.);
279   fLineExpTimeMax->SetLineColor(kGreen);
280   fLineExpTimeMax->SetLineWidth(2);
281   
282   fLineExpTotMin = new TLine( 5., 0., 5., 0.);
283   fLineExpTotMin->SetLineColor(kGreen);
284   fLineExpTotMin->SetLineWidth(2);
285   
286   fLineExpTotMax = new TLine(20., 0., 20., 0.);
287   fLineExpTotMax->SetLineColor(kGreen);
288   fLineExpTotMax->SetLineWidth(2);
289   
290   h5->GetListOfFunctions()->Add(fLineExpTimeMin);
291   h5->GetListOfFunctions()->Add(fLineExpTimeMax);
292   h10->GetListOfFunctions()->Add(fLineExpTotMin);
293   h10->GetListOfFunctions()->Add(fLineExpTotMax);
294   
295   for (Int_t sm=0;sm<10;sm++){
296       fLineSMid035[sm] = new TLine( 40*sm, 0, 40*sm, 0.);
297       fLineSMid035[sm]->SetLineColor(kMagenta);
298       fLineSMid035[sm]->SetLineWidth(2);
299       GetRawsData(16)->GetListOfFunctions()->Add(fLineSMid035[sm]);
300       fLineSMid3671[sm] = new TLine( 40*sm+360, 0, 40*sm+360, 0.);
301       fLineSMid3671[sm]->SetLineColor(kMagenta);
302       fLineSMid3671[sm]->SetLineWidth(2);
303       GetRawsData(17)->GetListOfFunctions()->Add(fLineSMid3671[sm]);
304   }
305   
306 }
307
308 //____________________________________________________________________________ 
309 void AliTOFQADataMakerRec::InitRecPoints()
310 {
311   //
312   // create RecPoints histograms in RecPoints subdir
313   //
314
315   const Bool_t expert   = kTRUE ; 
316   const Bool_t image    = kTRUE ; 
317
318   TH1I * h0 = new TH1I("hTOFRecPoints",    "TOF RecPoints multiplicity ; TOF RecPoints number;Counts",200, 0, 200) ;
319
320   TH1F * h1 = new TH1F("hTOFRecPointsTimeIA", "RecPoints Time Spectrum in TOF (ns)- I/A side; Calibrated TOF time [ns];Counts", 25000, 0., 610.) ; 
321   TH1F * h2 = new TH1F("hTOFRecPointsTimeOA", "RecPoints Time Spectrum in TOF (ns)- O/A side; Calibrated TOF time [ns];Counts", 25000, 0., 610.) ;
322   TH1F * h3 = new TH1F("hTOFRecPointsTimeIC", "RecPoints Time Spectrum in TOF (ns)- I/C side; Calibrated TOF time [ns];Counts", 25000, 0., 610.) ;
323   TH1F * h4 = new TH1F("hTOFRecPointsTimeOC", "RecPoints Time Spectrum in TOF (ns)- O/C side; Calibrated TOF time [ns];Counts", 25000, 0., 610.) ;
324   
325   TH1F * h5  = new TH1F("hTOFRecPointsRawTimeIA", "RecPoints raw Time Spectrum in TOF (ns)-I/A side; Measured TOF time [ns];Counts", 25000, 0., 610.) ; 
326   TH1F * h6  = new TH1F("hTOFRecPointsRawTimeOA", "RecPoints raw Time Spectrum in TOF (ns)-O/A side; Measured TOF time [ns];Counts", 25000, 0., 610.) ; 
327   TH1F * h7  = new TH1F("hTOFRecPointsRawTimeIC", "RecPoints raw Time Spectrum in TOF (ns)-I/C side; Measured TOF time [ns];Counts", 25000, 0., 610.) ; 
328   TH1F * h8  = new TH1F("hTOFRecPointsRawTimeOC", "RecPoints raw Time Spectrum in TOF (ns)-O/C side; Measured TOF time [ns];Counts", 25000, 0., 610.) ; 
329  
330   TH1F * h9   = new TH1F("hTOFRecPointsToTIA", "RecPoints ToT Spectrum in TOF (ns)-I/A side; Measured TOT [ns];Counts", 1000, 0., 48.8 ) ; 
331   TH1F * h10  = new TH1F("hTOFRecPointsToTOA", "RecPoints ToT Spectrum in TOF (ns)-O/A side; Measured TOT [ns];Counts", 1000, 0., 48.8 ) ; 
332   TH1F * h11  = new TH1F("hTOFRecPointsToTIC", "RecPoints ToT Spectrum in TOF (ns)-I/C side; Measured TOT [ns];Counts", 1000, 0., 48.8 ) ; 
333   TH1F * h12  = new TH1F("hTOFRecPointsToTOC", "RecPoints ToT Spectrum in TOF (ns)-O/C side; Measured TOT [ns];Counts", 1000, 0., 48.8 ) ; 
334   TH1F * h13 = new TH1F("hTOFRawChannelHits","TOF channel hits count; Channel ID; Hits", 158000, 0., 158000);
335
336   TH2F * h14 = new TH2F("hTOFRecPointsClusMap","RecPoints vs TOF phi-eta; eta (2*strip+padz);phi (48*sector+padx)",183, -0.5, 182.5,865,-0.5,864.5) ; 
337   TH2F * h15 = new TH2F("hTOFRecTimeVsStrip","RecPoints TOF time vs. strip (theta); Strip index; RecPoints TOF time (ns) ",92,-1.,91, 250, 0., 610.) ;
338
339   TH2F * h16 = new TH2F("hTOFRecPointsTimeVsTRM035","TOF RecPoints time vs TRM - crates 0 to 35; TRM index = DDL*10+TRM(0-9);TOF time [ns]", 361, 0., 361., 250, 0., 610.0) ;
340   TH2F * h17 = new TH2F("hTOFRecPointsTimeVsTRM3671","TOF RecPoints time vs TRM - crates 36 to 72; TRM index = DDL**10+TRM(0-9);TOF time [ns]", 361, 360., 721., 250, 0., 610.0) ;
341   TH2F * h18 = new TH2F("hTOFRecPointsToTVsTRM035","TOF RecPoints ToT vs TRM - crates 0 to 35; TRM index = DDL*10+TRM(0-9);TOF ToT [ns] ", 361, 0., 361., 100, 0., 48.8) ;
342   TH2F * h19 = new TH2F("hTOFRecPointsToTVsTRM3671","TOF RecPoints ToT vs TRM - crates 36 to 72; TRM index = DDL*10+TRM(0-9);TOF ToT [ns] ", 361, 360., 721., 100, 0., 48.8) ;
343  
344     h0->Sumw2() ;
345     h1->Sumw2() ;
346     h2->Sumw2() ;
347     h3->Sumw2() ;
348     h4->Sumw2() ;
349     h5->Sumw2() ;
350     h6->Sumw2() ;
351     h7->Sumw2() ;
352     h8->Sumw2() ;
353     h9->Sumw2() ;
354     h10->Sumw2() ;
355     h11->Sumw2() ;
356     h12->Sumw2() ;
357     h13->Sumw2() ;
358     h14->Sumw2() ;
359     h15->Sumw2() ;
360     h16->Sumw2() ;
361     h17->Sumw2() ;
362     h18->Sumw2() ;
363     h19->Sumw2() ;
364
365   Add2RecPointsList(h0, 0,   !expert,  image) ;
366   Add2RecPointsList(h1, 1,   !expert,  image) ;
367   Add2RecPointsList(h2, 2,   !expert,  image) ;
368   Add2RecPointsList(h3, 3,   !expert,  image) ;
369   Add2RecPointsList(h4, 4,   !expert,  image) ;
370   Add2RecPointsList(h5, 5,    expert, !image) ;
371   Add2RecPointsList(h6, 6,    expert, !image) ;
372   Add2RecPointsList(h7, 7,    expert, !image) ;
373   Add2RecPointsList(h8, 8,    expert, !image) ;
374   Add2RecPointsList(h9, 9,   !expert,  image) ;
375   Add2RecPointsList(h10, 10, !expert,  image) ;
376   Add2RecPointsList(h11, 11, !expert,  image) ;
377   Add2RecPointsList(h12, 12, !expert,  image) ;
378   Add2RecPointsList(h13, 13,  expert, !image) ;
379   Add2RecPointsList(h14, 14,  expert, !image) ;
380   Add2RecPointsList(h15, 15,  expert, !image) ;
381   Add2RecPointsList(h16, 16,  expert, !image) ;
382   Add2RecPointsList(h17, 17,  expert, !image) ;
383   Add2RecPointsList(h18, 18,  expert, !image) ;
384   Add2RecPointsList(h19, 19,  expert, !image) ;
385 }
386
387 //____________________________________________________________________________ 
388 void AliTOFQADataMakerRec::InitESDs()
389 {
390   //
391   //create ESDs histograms in ESDs subdir
392   //
393
394   const Bool_t expert   = kTRUE ; 
395   const Bool_t image    = kTRUE ; 
396
397   TH1I * h0  = new TH1I("hTOFESDs", "Number of matched TOF tracks per event;Number of TOF matched ESD tracks;Counts", 200, 0, 200) ;  
398   TH1F * h1  = new TH1F("hTOFESDsTime", "Matched  ESDs tracks: TOF Time spectrum; Calibrated TOF time [ns];Counts", 25000, 0., 610. ) ; 
399   TH1F * h2  = new TH1F("hTOFESDsRawTime", "Matched ESDs tracks: TOF raw Time spectrum;Measured TOF time [ns];Counts", 25000, 0., 610.) ; 
400   TH1F * h3  = new TH1F("hTOFESDsToT", "Matched ESDs tracks: TOF ToT spectrum; ESDs ToT [ns];Counts",1000, 0., 48.8) ; 
401   TH1F * h4  = new TH1F("hTOFESDsPIDoverMatched", "Fraction of TOF identified tracks over matched TOF tracks per event; identified by TOF / matched tracks [%];Counts", 101, -1., 100.) ;  
402   TH1F * h5  = new TH1F("hTOFESDsPID", "Fraction of TOF identified tracks over ESD identified tracks per event;  identified by TOF / ESD identified [%];Counts", 101, -1., 100.) ;  
403   TH1F * h6  = new TH1F("hTOFESDsTPCmatched", "TPC-TOF matched tracks momentum distribution (GeV/c); p (GeV/c);Counts", 50, 0.20, 5.00) ;  
404   TH1F * h7  = new TH1F("hTOFESDsMatchingProb", "TPC-TOF track-matching probability per event;TPC-TOF track-matching probability (%)  ;Counts",101, -1.0, 100) ;  
405   TH1F * h8  = new TH1F("hTOFESDsHitOverTracked", "Fraction of TOF matching tracks over propagated tracks per event; TOF matched / propagated to TOF [%];Counts",101, -1.0, 100.0) ;  
406   TH1F * h9  = new TH1F("hTOFESDsDiffTime", "ESDs t_{TOF}-t_{exp} spectrum in TOF (ps); t_{TOF}-t_{exp} [ps];Counts", 200, -2440., 2440.) ; 
407   TH1F * h10  = new TH1F("hTOFHitsLength", "Matched ESDs tracks: Length Spectrum; Track length [cm];Counts", 800, 0., 800) ; 
408   h0->Sumw2() ;
409   h1->Sumw2() ;
410   h2->Sumw2() ;
411   h3->Sumw2() ;
412   h4->Sumw2() ;
413   h5->Sumw2() ;
414   h6->Sumw2() ;
415   h7->Sumw2() ;
416   h8->Sumw2() ;
417   h9->Sumw2() ;
418   h10->Sumw2() ;
419
420   Add2ESDsList(h0, 0, !expert,  image) ;
421   Add2ESDsList(h1, 1, !expert,  image) ;
422   Add2ESDsList(h2, 2,  expert,  image) ;
423   Add2ESDsList(h3, 3, !expert,  image) ;
424   Add2ESDsList(h4, 4,  expert,  image) ;
425   Add2ESDsList(h5, 5,  expert,  image) ;
426   Add2ESDsList(h6, 6,  expert,  image) ; 
427   Add2ESDsList(h7, 7,  expert,  image) ; 
428   Add2ESDsList(h8, 8,  expert,  image) ; 
429   Add2ESDsList(h9, 9, !expert,  image) ;
430   Add2ESDsList(h10, 10, !expert,  image) ;
431 }
432
433
434 //____________________________________________________________________________
435 void AliTOFQADataMakerRec::MakeRaws(AliRawReader* rawReader)
436 {
437   //
438   // makes data from Raws
439   //
440     if (rawReader->GetType()==7) {
441         
442         Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;//in ns
443         Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
444         
445         Int_t ntof[5]; /* 0=tot, 1=IA, 2=OA, 3=IC, 4=OC*/
446         for (Int_t j=0;j<5;j++){ ntof[j]=0;}
447         
448         Int_t equipmentID[5]; //(ddl, trm, chain,tdc,channel)
449         Int_t volumeID[5];   //(sector,plate,strip,padX,padZ)
450         Int_t volumeID2[5];   //(sector,plate,strip,padZ,padX) to use AliTOFGeometry::GetIndex()
451         Int_t out[5]; //   out=(indexZ,indexPhi)   
452         Int_t chIndex=-1;
453         
454         TClonesArray * clonesRawData;
455         //AliTOFRawStream tofInput(rawReader);
456         fTOFRawStream.SetRawReader(rawReader);
457
458         //uncomment if needed to apply DeltaBC correction
459         //fTOFRawStream.ApplyBCCorrections(kTRUE);
460         
461         for (Int_t iDDL = 0; iDDL < AliTOFGeometry::NDDL()*AliTOFGeometry::NSectors(); iDDL++){
462             rawReader->Reset();
463             
464             fTOFRawStream.LoadRawDataBuffersV2(iDDL);
465             clonesRawData = (TClonesArray*)fTOFRawStream.GetRawData();
466             for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
467                 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
468                 
469                 if (tofRawDatum->GetTOF()){
470                     
471                     equipmentID[0]=iDDL;
472                     equipmentID[1]=tofRawDatum->GetTRM(); 
473                     equipmentID[2]=tofRawDatum->GetTRMchain();
474                     equipmentID[3]=tofRawDatum->GetTDC();
475                     equipmentID[4]=tofRawDatum->GetTDCchannel();
476                     
477                     if (CheckEquipID(equipmentID)){
478                         fTOFRawStream.EquipmentId2VolumeId(iDDL, 
479                                                            tofRawDatum->GetTRM(), 
480                                                            tofRawDatum->GetTRMchain(),
481                                                            tofRawDatum->GetTDC(), 
482                                                            tofRawDatum->GetTDCchannel(), 
483                                                            volumeID);
484                         if (FilterSpare(equipmentID)) continue;
485                         if (FilterLTMData(equipmentID)){ //counts LTM hits
486                             if (tofRawDatum->GetTOT()) GetRawsData(15)->Fill(equipmentID[0]);
487                         } else {
488                             if (CheckVolumeID(volumeID)){  
489                                 
490                                 GetMapIndeces(volumeID,out);
491                                 volumeID2[0]=volumeID[0];
492                                 volumeID2[1]=volumeID[1];
493                                 volumeID2[2]=volumeID[2];
494                                 volumeID2[3]=volumeID[4];
495                                 volumeID2[4]=volumeID[3];
496                                 chIndex=AliTOFGeometry::GetIndex(volumeID2);
497                                 
498                                 if (tofRawDatum->GetTOT()){         
499                                     if (!(fCalibData->GetNoiseStatus(chIndex)==AliTOFChannelOnlineStatusArray::kTOFNoiseBad)
500                                         && (fCalibData->GetHWStatus(chIndex) == AliTOFChannelOnlineStatusArray::kTOFHWOk)) {//noise and enabled filter
501                                         ntof[0]++; //counter for tof hits
502                                         
503                                         //fill global spectra for DQM plots
504                                         GetRawsData(5)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;//in ns
505                                         GetRawsData(10)->Fill( tofRawDatum->GetTOT()*tot2ns) ;//in ns
506                                         
507                                         //fill side-related spectra for experts plots
508                                         if (volumeID2[0]>4 && volumeID2[0]<14){       //I side
509                                             if ((iDDL%4==0)|| (iDDL%4==1)){ //A side
510                                                 ntof[1]++;
511                                                 GetRawsData(6)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;
512                                                 GetRawsData(11)->Fill( tofRawDatum->GetTOT()*tot2ns) ;
513                                             } else {
514                                                 if ((iDDL%4==2)|| (iDDL%4==3)){//C side
515                                                     ntof[3]++;
516                                                     GetRawsData(8)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;
517                                                     GetRawsData(13)->Fill( tofRawDatum->GetTOT()*tot2ns) ;
518                                                 }
519                                             }
520                                         } else {                                    
521                                             if (volumeID2[0]<5 || volumeID2[0]>13){   //O side
522                                                 if ((iDDL%4==0)|| (iDDL%4==1)){ //A side
523                                                     ntof[2]++;
524                                                     GetRawsData(7)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;
525                                                     GetRawsData(12)->Fill( tofRawDatum->GetTOT()*tot2ns) ;
526                                                 } else {
527                                                     if ((iDDL%4==2)|| (iDDL%4==3)){//C side
528                                                         ntof[4]++;
529                                                         GetRawsData(9)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;
530                                                         GetRawsData(14)->Fill( tofRawDatum->GetTOT()*tot2ns) ;
531                                                     }
532                                                 }
533                                             }   
534                                         }
535                                         GetRawsData(18)->Fill(chIndex);
536                                         //compute TRM offset
537                                         Int_t trm= iDDL*10+(equipmentID[1]-3);
538                                         if (iDDL>=0 && iDDL<36) {
539                                             GetRawsData(16)->Fill(trm) ;
540                                             GetRawsData(20)->Fill(trm,tofRawDatum->GetTOF()*tdc2ns);
541                                             GetRawsData(22)->Fill(trm,tofRawDatum->GetTOT()*tot2ns);
542                                         }
543                                         if (iDDL>=36 && iDDL<72) {
544                                             GetRawsData(17)->Fill(trm) ;//in ns 
545                                             GetRawsData(21)->Fill(trm,tofRawDatum->GetTOF()*tdc2ns);
546                                             GetRawsData(23)->Fill(trm,tofRawDatum->GetTOT()*tot2ns);
547                                         }                               
548                                         GetRawsData(24)->Fill(GetStripIndex(volumeID),tofRawDatum->GetTOF()*tdc2ns) ;
549                                         //GetRawsData(25)->Fill( out[0],out[1]) ;//raw map
550                                     }//noise filter
551                                 }//end hit selection
552                                 else { //orphans
553                                     if (!(fCalibData->GetNoiseStatus(chIndex) == AliTOFChannelOnlineStatusArray::kTOFNoiseBad)
554                                         && (fCalibData->GetHWStatus(chIndex) == AliTOFChannelOnlineStatusArray::kTOFHWOk))
555                                         GetRawsData(19)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;//in ns
556                                 }//end orphans
557                             }//end volumeID check
558                         }//end LTM filter
559                     }//end equipID check
560                 }
561             } //while loop
562             clonesRawData->Clear();
563         } // DDL Loop
564         
565         for (Int_t j=0;j<5;j++){
566           GetRawsData(j)->Fill(ntof[j]);
567         }
568         fProcessedRawEventN++;
569
570         fTOFRawStream.Clear();
571
572     } else {
573         AliDebug(1,Form("Event of type %d found. Skipping non-physics event for QA.\n", rawReader->GetType())); 
574     }
575     EnableDqmShifterOpt(kTRUE);
576 }
577
578 //____________________________________________________________________________
579 void AliTOFQADataMakerRec::MakeRecPoints(TTree * clustersTree)
580 {
581     //
582   // Make data from Clusters
583   //
584  
585   Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
586   Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
587
588   Int_t volumeID[5];
589   Int_t volumeID2[5];
590   Int_t out[5];
591
592   TBranch *branch=clustersTree->GetBranch("TOF");
593   if (!branch) { 
594     AliError("can't get the branch with the TOF clusters !");
595     return;
596   }
597
598   static TClonesArray dummy("AliTOFcluster",10000);
599   dummy.Clear();
600   TClonesArray *clusters=&dummy;
601   branch->SetAddress(&clusters);
602
603   // Import the tree
604   clustersTree->GetEvent(0);  
605  
606   GetRecPointsData(0)->Fill((Int_t)clusters->GetEntriesFast()) ; 
607   
608   TIter next(clusters) ; 
609   AliTOFcluster * c ; 
610   while ( (c = dynamic_cast<AliTOFcluster *>(next())) ) {
611    
612       volumeID[0] = c->GetDetInd(0);
613       volumeID[1] = c->GetDetInd(1);
614       volumeID[2] = c->GetDetInd(2);
615       volumeID[3] = c->GetDetInd(4); //X and Z indeces inverted in RecPoints
616       volumeID[4] = c->GetDetInd(3); //X and Z indeces inverted in RecPoints
617       
618       for (Int_t i=0;i<5;i++){
619           volumeID2[i]=c->GetDetInd(i); //X and Z indeces inverted in RecPoints
620       }
621       
622       GetMapIndeces(volumeID,out);
623       Int_t chIndex=AliTOFGeometry::GetIndex(volumeID2);
624       Int_t iDDL=AliTOFRawStream::Geant2DDL(volumeID2);
625       Int_t iTRM=AliTOFRawStream::Geant2TRM(volumeID2);
626
627       if (c->GetTDCRAW() && c->GetTDC() && c->GetToT()){
628           if (volumeID2[0]>4 && volumeID2[0]<14){       //I side
629               if ((iDDL%4==0)|| (iDDL%4==1)){ //A side
630                   GetRecPointsData(1)->Fill( c->GetTDC()*tdc2ns) ;//in ns
631                   GetRecPointsData(5)->Fill( c->GetTDCRAW()*tdc2ns) ;//in ns
632                   GetRecPointsData(9)->Fill( c->GetToT()*tot2ns) ;//in ns
633                   
634               } else {
635                   if ((iDDL%4==2)|| (iDDL%4==3)){//C side
636                       GetRecPointsData(3)->Fill( c->GetTDC()*tdc2ns) ;//in ns
637                       GetRecPointsData(7)->Fill( c->GetTDCRAW()*tdc2ns) ;//in ns
638                       GetRecPointsData(11)->Fill( c->GetToT()*tot2ns) ;//in ns
639                   }
640               }
641           } else {
642               if (volumeID2[0]<5 || volumeID2[0]>13){       //O side
643                   if ((iDDL%4==0)|| (iDDL%4==1)){ //A side
644                       GetRecPointsData(2)->Fill( c->GetTDC()*tdc2ns) ;//in ns
645                       GetRecPointsData(6)->Fill( c->GetTDCRAW()*tdc2ns) ;//in ns
646                       GetRecPointsData(10)->Fill( c->GetToT()*tot2ns) ;//in ns
647                   } else {
648                       if ((iDDL%4==2)|| (iDDL%4==3)){//C side
649                           GetRecPointsData(4)->Fill( c->GetTDC()*tdc2ns) ;//in ns
650                           GetRecPointsData(8)->Fill( c->GetTDCRAW()*tdc2ns) ;//in ns
651                           GetRecPointsData(12)->Fill( c->GetToT()*tot2ns) ;//in ns
652                       }
653                   }
654               }
655           }
656           GetRecPointsData(14)->Fill(out[0],out[1]);
657           GetRecPointsData(15)->Fill(GetStripIndex(volumeID), c->GetTDC()*tdc2ns) ;
658           Int_t trm= iDDL*10+(iTRM-3);
659           if (iDDL>=0 && iDDL<36) {
660               GetRecPointsData(16)->Fill(trm,c->GetTDC()*tdc2ns);
661               GetRecPointsData(18)->Fill(trm,c->GetToT()*tot2ns);
662           }
663           if (iDDL>=36 && iDDL<72) {
664               GetRecPointsData(17)->Fill(trm,c->GetTDC()*tdc2ns);
665               GetRecPointsData(19)->Fill(trm,c->GetToT()*tot2ns);
666           }
667           GetRecPointsData(13)->Fill(chIndex) ;//in ns
668       }//hit selection
669   }//end while  
670  
671   EnableDqmShifterOpt(kFALSE);
672 }
673
674 //____________________________________________________________________________
675 void AliTOFQADataMakerRec::MakeESDs(AliESDEvent * esd)
676 {
677   //
678   // make QA data from ESDs
679   //  
680     const Double_t speedOfLight = TMath::C()*1E2*1E-12; // cm/ps
681     const Double_t pionMass = 0.13957018; //GeV/c^2
682
683     Int_t ntrk = esd->GetNumberOfTracks() ; 
684     Int_t ntof=0;
685     Int_t ntofpid=0;
686     Int_t nesdpid=0;
687     Int_t ntpc=0;
688     Int_t ntpctof=0;
689     Int_t ntofout=0;
690     
691     while (ntrk--) {
692         AliESDtrack *track=esd->GetTrack(ntrk);
693         Double_t tofTime=track->GetTOFsignal();//in ps
694         Double_t tofTimeRaw=track->GetTOFsignalRaw();//in ps
695         Double_t tofToT=track->GetTOFsignalToT(); //in ps
696
697         UInt_t status=track->GetStatus();
698         if (track->IsOn(AliESDtrack::kTPCrefit)) {
699             ntpc++;
700             Double_t y=track->Eta();
701             if ((status&AliESDtrack::kTOFout)!=0) {
702                 ntofout++;
703                 if (tofTime>0){
704                     ntof++;
705                     
706                     if (TMath::Abs(y)<0.9) {
707                         GetESDsData(6)->Fill(track->Pt());
708                         ntpctof++;
709                     }
710                     GetESDsData(1)->Fill(tofTime*1E-3);
711                     GetESDsData(2)->Fill(tofTimeRaw*1E-3); 
712                     GetESDsData(3)->Fill(tofToT*1E-3);
713                     //check how many tracks where ESD PID is ok 
714                     if ((status&AliESDtrack::kESDpid)!=0) nesdpid++; 
715                     if (((status&AliESDtrack::kESDpid)&AliESDtrack::kTOFpid)!=0) ntofpid++;
716                     
717                     Double_t length =track->GetIntegratedLength();
718                     Double_t mom2=(track->Pt()*track->Pt())+(track->Pz()*track->Pz());
719                     Double_t eTexp = TMath::Sqrt(1+(pionMass*pionMass/mom2))*length/speedOfLight; //in ps
720                     GetESDsData(9)->Fill(tofTime-eTexp);
721                     GetESDsData(10)->Fill(length);
722                 } //end check on matched tracks
723             }
724         }//end check on TPCrefit
725     }
726     
727     GetESDsData(0)->Fill(ntof) ;
728   
729     
730     if(ntof>0) {
731         Float_t ratio = (Int_t)ntofpid/(Int_t)ntof*100; //identified by TOF over matched
732         GetESDsData(4)->Fill(ratio) ;
733     }
734     
735     if(nesdpid>0) {
736         Float_t ratio = (Float_t)ntofpid/(Float_t)nesdpid *100; //identified by TOF over identified
737         GetESDsData(5)->Fill(ratio) ;
738     }
739     
740     if(ntpc>0){
741         Float_t ratio = (Float_t)ntof/(Float_t)ntpc*100.; //matching probability
742         GetESDsData(7)->Fill(ratio) ;
743     }
744     
745     if(ntofout>0) {
746         Float_t ratio = (Float_t)ntof/(Float_t)ntofout*100; //matched over propagated to TOF outer radius
747         GetESDsData(8)->Fill(ratio) ;
748     }
749     EnableDqmShifterOpt(kFALSE);
750 }
751
752 //____________________________________________________________________________ 
753 void AliTOFQADataMakerRec::StartOfDetectorCycle()
754 {
755   //
756   //Detector specific actions at start of cycle
757  
758   fCalibData = GetCalibData();
759
760 }
761
762 //____________________________________________________________________________ 
763 void AliTOFQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
764 {
765   //Detector specific actions at end of cycle
766   // do the QA checking
767     for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
768         if ( !AliQAv1::Instance()->IsEventSpecieSet(specie) ) 
769             continue ; 
770         
771         if (fEnableDqmShifterOpt){
772             // Help make the raw qa histogram easier to interpret for the DQM shifter
773             if (!GetRawsData(0) || !GetRawsData(5) || !GetRawsData(10) 
774                 || !GetRawsData(15) || !GetRawsData(16) || !GetRawsData(17)) {
775                 printf("No histogram for DQM found - Possible memory corruption ???. Please check\n") ; 
776                 continue;
777             }
778             AliInfo(Form("Processed %i physics raw events",fProcessedRawEventN));
779             
780             
781             if (fCalibData){
782                 //normalize TRM hits plots to the number of enabled channels from OCDB object
783                 TH1F * hTrmChannels035 = new TH1F("hTrmchannels035", "Active channels per TRM - crates 0 to 35;TRM index = SMid(crate*10)+TRM(0-9);Active channels",  361, 0., 361.) ;
784                 TH1F * hTrmChannels3671 = new TH1F("hTrmChannels3671","Active channels per TRM - crates 36 to 71 ;TRM index = SMid(crate*10)+TRM(0-9);Active channels", 361, 360., 721.) ;
785                 
786                 for (Int_t ch = 0; ch <  fCalibData->GetSize(); ch++) {
787                     if (!(fCalibData->GetNoiseStatus(ch)==AliTOFChannelOnlineStatusArray::kTOFNoiseBad)
788                         && (fCalibData->GetHWStatus(ch) == AliTOFChannelOnlineStatusArray::kTOFHWOk)){
789                         Int_t geoId[5];
790                         Int_t detId[5];
791                         AliTOFGeometry::GetVolumeIndices(ch,geoId);
792                         AliTOFRawStream::Geant2EquipmentId(geoId,detId); //detID=(ddl,trm,tdc, chain,channel)
793                         if (detId[0]<36) hTrmChannels035->Fill((detId[1]-3)+detId[0]*10);
794                         else hTrmChannels3671->Fill((detId[1]-3)+detId[0]*10);
795                     }
796                 }
797                 GetRawsData(16)->Divide(hTrmChannels035);
798                 GetRawsData(16)->SetTitle("TRMs average hit number per active channel - crates 0-35");
799                 GetRawsData(16)->GetYaxis()->SetTitle("hits/active channels");
800                 GetRawsData(17)->Divide(hTrmChannels3671);
801                 GetRawsData(17)->SetTitle("TRMs average hit number per active channel - crates 36-71");
802                 GetRawsData(17)->GetYaxis()->SetTitle("hits/active channels");
803             }
804             
805             //set minima and maxima to allow log scale
806             Double_t yTimeMax = GetRawsData(5)->GetMaximum();
807             Double_t yTotMax = GetRawsData(10)->GetMaximum();
808             Double_t yTrmMax = TMath::Max(GetRawsData(16)->GetMaximum(),GetRawsData(17)->GetMaximum());
809             // Double_t yLtmMax = GetRawsData(15)->GetMaximum();
810             // Double_t yHitMax = GetRawsData(0)->GetMaximum();
811
812             // GetRawsData(0)->SetMinimum(0.1);
813             // GetRawsData(5)->SetMinimum(0.1);
814             // GetRawsData(10)->SetMinimum(0.1);
815             // GetRawsData(15)->SetMinimum(0.1);
816             GetRawsData(16)->SetMinimum(0.0001);
817             GetRawsData(17)->SetMinimum(0.0001);
818
819             // GetRawsData(0)->SetMaximum(yHitMax);
820             // GetRawsData(5)->SetMaximum(yTimeMax);
821             // GetRawsData(10)->SetMaximum(yTotMax);
822             // GetRawsData(15)->SetMaximum(yLtmMax);
823             GetRawsData(16)->SetMaximum(yTrmMax);
824             GetRawsData(17)->SetMaximum(yTrmMax);
825             
826             fLineExpTimeMin->SetY2(yTimeMax);
827             fLineExpTimeMax->SetY2(yTimeMax);
828             fLineExpTotMin->SetY2(yTotMax);
829             fLineExpTotMax->SetY2(yTotMax);
830             for (Int_t sm=0;sm<10;sm++){
831                 fLineSMid035[sm]->SetY2(yTrmMax);
832                 fLineSMid3671[sm]->SetY2(yTrmMax);
833             }
834
835
836             //make up for all histos 
837             for(Int_t j=0;j<20;j++){
838                 if (j<5) {
839                     GetRawsData(j)->SetMarkerColor(kRed);
840                     GetRawsData(j)->SetMarkerStyle(7);
841                 } else {
842                     GetRawsData(j)->SetLineColor(kBlue+1);
843                     GetRawsData(j)->SetLineWidth(1);
844                     GetRawsData(j)->SetMarkerColor(kBlue+1);
845                 }
846             }
847             
848             for (Int_t j=15;j<19;j++){
849                 GetRawsData(j)->SetFillColor(kGray+1);
850                 GetRawsData(j)->SetOption("bar");
851             }
852             
853         }
854     }
855     AliQAChecker::Instance()->Run(AliQAv1::kTOF, task, list) ;  
856 }
857 //____________________________________________________________________________
858 void AliTOFQADataMakerRec::GetMapIndeces(const Int_t* const in , Int_t* out)
859 {
860   //
861   //return appropriate indeces for the theta-phi map
862   //
863
864   Int_t npadX = AliTOFGeometry::NpadX();
865   Int_t npadZ = AliTOFGeometry::NpadZ();
866   Int_t nStripA = AliTOFGeometry::NStripA();
867   Int_t nStripB = AliTOFGeometry::NStripB();
868   Int_t nStripC = AliTOFGeometry::NStripC();
869
870   Int_t isector = in[0];
871   Int_t iplate = in[1];
872   Int_t istrip = in[2];
873   Int_t ipadX = in[3]; 
874   Int_t ipadZ = in[4]; 
875   
876   Int_t stripOffset = 0;
877   switch (iplate) {
878   case 0:
879     stripOffset = 0;
880       break;
881   case 1:
882     stripOffset = nStripC;
883     break;
884   case 2:
885     stripOffset = nStripC+nStripB;
886     break;
887   case 3:
888     stripOffset = nStripC+nStripB+nStripA;
889     break;
890   case 4:
891     stripOffset = nStripC+nStripB+nStripA+nStripB;
892     break;
893   default:
894     AliDebug(1,Form("Wrong plate number in TOF (%d) !",iplate));
895     break;
896   };
897   Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
898   Int_t phiindex=npadX*isector+ipadX+1;
899   out[0]=zindex;  
900   out[1]=phiindex;  
901   
902 }
903
904 //---------------------------------------------------------------
905 Int_t AliTOFQADataMakerRec::GetStripIndex(const Int_t * const in)
906 {
907     /* return tof strip index between 0 and 91 */
908
909   Int_t nStripA = AliTOFGeometry::NStripA();
910   Int_t nStripB = AliTOFGeometry::NStripB();
911   Int_t nStripC = AliTOFGeometry::NStripC();
912
913   // Int_t isector = in[0];
914   Int_t iplate = in[1];
915   Int_t istrip = in[2];
916   //Int_t ipadX = in[3]; 
917   //Int_t ipadZ = in[4]; 
918   
919   Int_t stripOffset = 0;
920   switch (iplate) {
921   case 0:
922     stripOffset = 0;
923       break;
924   case 1:
925     stripOffset = nStripC;
926     break;
927   case 2:
928     stripOffset = nStripC+nStripB;
929     break;
930   case 3:
931     stripOffset = nStripC+nStripB+nStripA;
932     break;
933   case 4:
934     stripOffset = nStripC+nStripB+nStripA+nStripB;
935     break;
936   default:
937       AliDebug(1,Form("Wrong plate number in TOF (%d) !",iplate));
938       stripOffset=-1;
939       break;
940   };
941   
942   if (stripOffset<0 || stripOffset>92) return -1;
943   else 
944       return (stripOffset+istrip);
945   
946 }
947 //---------------------------------------------------------------
948 Bool_t  AliTOFQADataMakerRec::CheckVolumeID(const Int_t * const volumeID)
949 {
950     //
951     //Checks volume ID validity
952     //
953    
954     for (Int_t j=0;j<5;j++){
955         if (volumeID[j]<0) {
956             AliDebug(1,Form("Invalid detector volume index for volumeID[%i]",j));
957             return kFALSE;
958         }
959     }
960     return kTRUE;
961     
962 }
963
964 //---------------------------------------------------------------
965 Bool_t  AliTOFQADataMakerRec::CheckEquipID(const Int_t * const equipmentID)
966 {
967     //
968     //Checks equipment ID validity
969     
970    for (Int_t j=0;j<5;j++){
971         if (equipmentID[j]<0) {
972           AliDebug(1,Form("Invalid equipment volume index for equipmentID[%i]",j));
973           return kFALSE;
974         }
975    }
976    return kTRUE;
977 }
978 //---------------------------------------------------------------
979 Bool_t  AliTOFQADataMakerRec::FilterLTMData(const Int_t * const equipmentID) const
980 {
981   /*It returns kTRUE if data come from LTM.
982     It thus filters trigger-related signals  */
983
984   Int_t ddl, trm, tdc;
985   //if (!CheckEquipID(equipmentID)) return kFALSE;
986   ddl = equipmentID[0];
987   trm = equipmentID[1];
988   tdc = equipmentID[3];
989   
990   if ((ddl%2==1) && (trm==3) && (tdc>11 && tdc<15))
991     return kTRUE;
992   else 
993     return kFALSE;
994  
995 }
996 //---------------------------------------------------------------
997 Bool_t  AliTOFQADataMakerRec::FilterSpare(const Int_t * const equipmentID) const
998 {
999   /*It returns kTRUE if data come from spare 
1000     equipment ID. 
1001     So far only check on TRM 3 crate left is implemented */
1002
1003   Int_t ddl, trm, tdc;
1004   //if (!CheckEquipID(equipmentID)) return kFALSE;
1005   ddl = equipmentID[0];
1006   trm = equipmentID[1];
1007   tdc = equipmentID[3];
1008   
1009   if ((ddl%2==1) && (trm==3) && (tdc>2 && tdc<12))
1010     return kTRUE;
1011   else 
1012     return kFALSE;
1013  
1014 }