]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFQADataMakerRec.cxx
bugfix: delete old deprecated duplication of code, custom component libraries had...
[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 /* Modified by fbellini on 22/04/2010
25    - Added filter for physics events
26
27    Modified by fbellini on 16/04/2010
28    - Added EnableDqmShifterOpt() 
29    - Modified EndOfDetectorCycle() with options for DQM         
30    - Updated ESDs QA
31
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
37
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  
43    - Updated Raw QA
44 */
45
46 #include <TClonesArray.h>
47 #include <TH1F.h> 
48 #include <TH2F.h> 
49 #include <TLine.h>
50
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"
63
64
65 ClassImp(AliTOFQADataMakerRec)
66            
67 //____________________________________________________________________________ 
68   AliTOFQADataMakerRec::AliTOFQADataMakerRec() : 
69   AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kTOF), "TOF Quality Assurance Data Maker"),
70   fCalibData(0x0),
71   fEnableNoiseFiltering(kFALSE),
72     fEnableDqmShifterOpt(kFALSE),
73     fProcessedRawEventN(0)
74 {
75   //
76   // ctor
77   //
78    
79 }
80
81 //____________________________________________________________________________ 
82 AliTOFQADataMakerRec::AliTOFQADataMakerRec(const AliTOFQADataMakerRec& qadm) :
83   AliQADataMakerRec(),
84   fCalibData(qadm.fCalibData),
85   fEnableNoiseFiltering(qadm.fEnableNoiseFiltering),
86   fEnableDqmShifterOpt(qadm.fEnableDqmShifterOpt),
87   fProcessedRawEventN(qadm.fProcessedRawEventN)
88 {
89   //
90   //copy ctor 
91   //
92   SetName((const char*)qadm.GetName()) ; 
93   SetTitle((const char*)qadm.GetTitle()); 
94   
95  
96 }
97
98 //__________________________________________________________________
99 AliTOFQADataMakerRec& AliTOFQADataMakerRec::operator = (const AliTOFQADataMakerRec& qadm )
100 {
101   //
102   // assignment operator.
103   //
104   this->~AliTOFQADataMakerRec();
105   new(this) AliTOFQADataMakerRec(qadm);
106   return *this;
107 }
108  
109 //----------------------------------------------------------------------------
110 AliTOFChannelOnlineStatusArray* AliTOFQADataMakerRec::GetCalibData() const
111 {
112   //
113   // Retrive TOF calib objects from OCDB
114   //
115     AliCDBManager *man = AliCDBManager::Instance();
116
117     AliCDBEntry *cdbe=0;
118
119     cdbe = man->Get("TOF/Calib/Status",fRun);
120     if(!cdbe){
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);
125     }
126     // Retrieval of data in directory TOF/Calib/Data:
127     
128     AliTOFChannelOnlineStatusArray * array = 0;
129     if (cdbe) array = (AliTOFChannelOnlineStatusArray *)cdbe->GetObject();
130     if (!array)  AliFatal("No calibration data from calibration database !");
131     
132     return array;
133     
134 }
135
136
137
138 //____________________________________________________________________________ 
139 void AliTOFQADataMakerRec::InitRaws()
140 {
141   //
142   // create Raws histograms in Raws subdir
143   //
144   const Bool_t expert   = kTRUE ; 
145   const Bool_t saveCorr = kTRUE ; 
146   const Bool_t image    = kTRUE ; 
147
148   TH1I * h0 =  new TH1I("hTOFRaws",      "TOF raw hit multiplicity; TOF raw hits number;Counts ",200, 0, 200) ;
149
150   TH1I * h1 =  new TH1I("hTOFRawsIA",    "TOF raw hit multiplicity - I/A side; TOF raw hits number ;Counts ",200, 0, 200) ;
151   TH1I * h2 =  new TH1I("hTOFRawsOA",    "TOF raw hit multiplicity - O/A side; TOF raw hits number ;Counts ",200, 0, 200) ;
152   TH1I * h3 =  new TH1I("hTOFRawsIC",    "TOF raw hit multiplicity - I/C side; TOF raw hits number ;Counts ",200, 0, 200) ;
153   TH1I * h4 =  new TH1I("hTOFRawsOC",    "TOF raw hit multiplicity - O/C side; TOF raw hits number ;Counts ",200, 0, 200) ;
154
155   TH1F * h5  = new TH1F("hTOFRawsTime", "TOF Raws - Hit time (ns);Measured Hit time [ns];Counts", 25000,0. ,610.) ; 
156
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.) ; 
161
162   TH1F * h10  = new TH1F("hTOFRawsToT", "TOF Raws - Hit ToT (ns);Measured Hit ToT (ns);Counts", 1000, 0., 48.8) ; 
163
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) ; 
168
169   TH1I * h15 = new TH1I("hTOFRawsLTMHits", "LTM hits ; Crate; Counts",  72, 0., 72.);
170   TH1I * h16  = new TH1I("hTOFRawsTRMHits035", "TRM hits  - crates 0 to 35 ;TRM index = SMid(crate*10)+TRM(0-9);Hits",  361, 0., 361.) ;
171   TH1I * h17  = new TH1I("hTOFRawsTRMHits3671","TRM hits  - crates 36 to 71 ;TRM index = SMid(crate*10)+TRM(0-9);Hits", 361, 360., 721.) ;
172   
173   TH1I * h18 = new TH1I("hTOFRawChannelHits","TOF channel hits count; Channel ID; Hits", 158000, 0., 158000);
174   
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.) ; 
181   
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) ; 
183
184   h0->Sumw2() ;
185   h1->Sumw2() ;
186   h2->Sumw2() ;
187   h3->Sumw2() ;
188   h4->Sumw2() ;
189   h5->Sumw2() ;
190   h6->Sumw2() ;
191   h7->Sumw2() ;
192   h8->Sumw2() ;
193   h9->Sumw2() ;
194   h10->Sumw2() ;
195   h11->Sumw2() ;
196   h12->Sumw2() ;
197   h13->Sumw2() ;
198   h14->Sumw2() ;
199   h15->Sumw2() ;
200   h16->Sumw2() ;
201   h17->Sumw2() ;
202   h18->Sumw2() ;
203   h19->Sumw2() ;
204   h20->Sumw2() ;
205   h21->Sumw2() ;
206   h22->Sumw2() ;
207   h23->Sumw2() ;
208   h24->Sumw2() ;
209
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) ;
235
236 }
237
238 //____________________________________________________________________________ 
239 void AliTOFQADataMakerRec::InitRecPoints()
240 {
241   //
242   // create RecPoints histograms in RecPoints subdir
243   //
244
245   const Bool_t expert   = kTRUE ; 
246   const Bool_t image    = kTRUE ; 
247
248   TH1I * h0 = new TH1I("hTOFRecPoints",    "TOF RecPoints multiplicity ; TOF RecPoints number;Counts",200, 0, 200) ;
249
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.) ;
254   
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.) ; 
259  
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);
265
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.) ;
268
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) ;
273  
274     h0->Sumw2() ;
275     h1->Sumw2() ;
276     h2->Sumw2() ;
277     h3->Sumw2() ;
278     h4->Sumw2() ;
279     h5->Sumw2() ;
280     h6->Sumw2() ;
281     h7->Sumw2() ;
282     h8->Sumw2() ;
283     h9->Sumw2() ;
284     h10->Sumw2() ;
285     h11->Sumw2() ;
286     h12->Sumw2() ;
287     h13->Sumw2() ;
288     h14->Sumw2() ;
289     h15->Sumw2() ;
290     h16->Sumw2() ;
291     h17->Sumw2() ;
292     h18->Sumw2() ;
293     h19->Sumw2() ;
294
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) ;
315 }
316
317 //____________________________________________________________________________ 
318 void AliTOFQADataMakerRec::InitESDs()
319 {
320   //
321   //create ESDs histograms in ESDs subdir
322   //
323
324   const Bool_t expert   = kTRUE ; 
325   const Bool_t image    = kTRUE ; 
326
327   TH1I * h0  = new TH1I("hTOFESDs", "Number of matched TOF tracks per event;Number of TOF matched ESD tracks;Counts", 200, 0, 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) ; 
338   h0->Sumw2() ;
339   h1->Sumw2() ;
340   h2->Sumw2() ;
341   h3->Sumw2() ;
342   h4->Sumw2() ;
343   h5->Sumw2() ;
344   h6->Sumw2() ;
345   h7->Sumw2() ;
346   h8->Sumw2() ;
347   h9->Sumw2() ;
348   h10->Sumw2() ;
349
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) ;
361 }
362
363
364 //____________________________________________________________________________
365 void AliTOFQADataMakerRec::MakeRaws(AliRawReader* rawReader)
366 {
367   //
368   // makes data from Raws
369   //
370     if (rawReader->GetType()==7) {
371         
372         Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;//in ns
373         Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
374         
375         Int_t ntof[5]; /* 0=tot, 1=IA, 2=OA, 3=IC, 4=OC*/
376         for (Int_t j=0;j<5;j++){ ntof[j]=0;}
377         
378         Int_t equipmentID[5]; //(ddl, trm, chain,tdc,channel)
379         Int_t volumeID[5];   //(sector,plate,strip,padX,padZ)
380         Int_t volumeID2[5];   //(sector,plate,strip,padZ,padX) to use AliTOFGeometry::GetIndex()
381         Int_t out[5]; //   out=(indexZ,indexPhi)   
382         Int_t chIndex=-1;
383         
384         TClonesArray * clonesRawData;
385         AliTOFRawStream tofInput(rawReader);
386         
387         //uncomment if needed to apply DeltaBC correction
388         //tofInput.ApplyBCCorrections(kTRUE);
389         
390         for (Int_t iDDL = 0; iDDL < AliTOFGeometry::NDDL()*AliTOFGeometry::NSectors(); iDDL++){
391             rawReader->Reset();
392             
393             tofInput.LoadRawDataBuffers(iDDL);
394             clonesRawData = (TClonesArray*)tofInput.GetRawData();
395             for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
396                 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
397                 
398                 if (tofRawDatum->GetTOF()){
399                     
400                     equipmentID[0]=iDDL;
401                     equipmentID[1]=tofRawDatum->GetTRM(); 
402                     equipmentID[2]=tofRawDatum->GetTRMchain();
403                     equipmentID[3]=tofRawDatum->GetTDC();
404                     equipmentID[4]=tofRawDatum->GetTDCchannel();
405                     
406                     if (CheckEquipID(equipmentID)){
407                         tofInput.EquipmentId2VolumeId(iDDL, 
408                                                       tofRawDatum->GetTRM(), 
409                                                       tofRawDatum->GetTRMchain(),
410                                                       tofRawDatum->GetTDC(), 
411                                                       tofRawDatum->GetTDCchannel(), 
412                                                       volumeID);
413                         if (FilterSpare(equipmentID)) continue;
414                         if (FilterLTMData(equipmentID)){ //counts LTM hits
415                             if (tofRawDatum->GetTOT()) GetRawsData(15)->Fill(equipmentID[0]);
416                         } else {
417                             if (CheckVolumeID(volumeID)){  
418                                 
419                                 GetMapIndeces(volumeID,out);
420                                 volumeID2[0]=volumeID[0];
421                                 volumeID2[1]=volumeID[1];
422                                 volumeID2[2]=volumeID[2];
423                                 volumeID2[3]=volumeID[4];
424                                 volumeID2[4]=volumeID[3];
425                                 chIndex=AliTOFGeometry::GetIndex(volumeID2);
426                                 
427                                 if (tofRawDatum->GetTOT()){         
428                                     if (!(fCalibData->GetNoiseStatus(chIndex)==AliTOFChannelOnlineStatusArray::kTOFNoiseBad)
429                                         && (fCalibData->GetHWStatus(chIndex) == AliTOFChannelOnlineStatusArray::kTOFHWOk)) {//noise and enabled filter
430                                         ntof[0]++; //counter for tof hits
431                                         
432                                         //fill global spectra for DQM plots
433                                         GetRawsData(5)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;//in ns
434                                         GetRawsData(10)->Fill( tofRawDatum->GetTOT()*tot2ns) ;//in ns
435                                         
436                                         //fill side-related spectra for experts plots
437                                         if (volumeID2[0]>4 && volumeID2[0]<14){       //I side
438                                             if ((iDDL%4==0)|| (iDDL%4==1)){ //A side
439                                                 ntof[1]++;
440                                                 GetRawsData(6)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;
441                                                 GetRawsData(11)->Fill( tofRawDatum->GetTOT()*tot2ns) ;
442                                             } else {
443                                                 if ((iDDL%4==2)|| (iDDL%4==3)){//C side
444                                                     ntof[3]++;
445                                                     GetRawsData(8)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;
446                                                     GetRawsData(13)->Fill( tofRawDatum->GetTOT()*tot2ns) ;
447                                                 }
448                                             }
449                                         } else {                                    
450                                             if (volumeID2[0]<5 || volumeID2[0]>13){   //O side
451                                                 if ((iDDL%4==0)|| (iDDL%4==1)){ //A side
452                                                     ntof[2]++;
453                                                     GetRawsData(7)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;
454                                                     GetRawsData(12)->Fill( tofRawDatum->GetTOT()*tot2ns) ;
455                                                 } else {
456                                                     if ((iDDL%4==2)|| (iDDL%4==3)){//C side
457                                                         ntof[4]++;
458                                                         GetRawsData(9)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;
459                                                         GetRawsData(14)->Fill( tofRawDatum->GetTOT()*tot2ns) ;
460                                                     }
461                                                 }
462                                             }   
463                                         }
464                                         GetRawsData(18)->Fill(chIndex);
465                                         //compute TRM offset
466                                         Int_t trm= iDDL*10+(equipmentID[1]-3);
467                                         if (iDDL>=0 && iDDL<36) {
468                                             GetRawsData(16)->Fill(trm) ;
469                                             GetRawsData(20)->Fill(trm,tofRawDatum->GetTOF()*tdc2ns);
470                                             GetRawsData(22)->Fill(trm,tofRawDatum->GetTOT()*tot2ns);
471                                         }
472                                         if (iDDL>=36 && iDDL<72) {
473                                             GetRawsData(17)->Fill(trm) ;//in ns 
474                                             GetRawsData(21)->Fill(trm,tofRawDatum->GetTOF()*tdc2ns);
475                                             GetRawsData(23)->Fill(trm,tofRawDatum->GetTOT()*tot2ns);
476                                         }                               
477                                         GetRawsData(24)->Fill(GetStripIndex(volumeID),tofRawDatum->GetTOF()*tdc2ns) ;
478                                         //GetRawsData(25)->Fill( out[0],out[1]) ;//raw map
479                                     }//noise filter
480                                 }//end hit selection
481                                 else { //orphans
482                                     if (!(fCalibData->GetNoiseStatus(chIndex) == AliTOFChannelOnlineStatusArray::kTOFNoiseBad)
483                                         && (fCalibData->GetHWStatus(chIndex) == AliTOFChannelOnlineStatusArray::kTOFHWOk))
484                                         GetRawsData(19)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;//in ns
485                                 }//end orphans
486                             }//end volumeID check
487                         }//end LTM filter
488                     }//end equipID check
489                 }
490             } //while loop
491             clonesRawData->Clear();
492         } // DDL Loop
493         
494         for (Int_t j=0;j<5;j++){
495             GetRawsData(j)->Fill(ntof[j]);
496         }
497         fProcessedRawEventN++;
498         
499     } else {
500         AliDebug(1,Form("Event of type %d found. Skipping non-physics event for QA.\n", rawReader->GetType())); 
501     }
502     EnableDqmShifterOpt(kTRUE);
503 }
504
505 //____________________________________________________________________________
506 void AliTOFQADataMakerRec::MakeRecPoints(TTree * clustersTree)
507 {
508     //
509   // Make data from Clusters
510   //
511  
512   Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
513   Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
514
515   Int_t volumeID[5];
516   Int_t volumeID2[5];
517   Int_t out[5];
518
519   TBranch *branch=clustersTree->GetBranch("TOF");
520   if (!branch) { 
521     AliError("can't get the branch with the TOF clusters !");
522     return;
523   }
524
525   static TClonesArray dummy("AliTOFcluster",10000);
526   dummy.Clear();
527   TClonesArray *clusters=&dummy;
528   branch->SetAddress(&clusters);
529
530   // Import the tree
531   clustersTree->GetEvent(0);  
532  
533   GetRecPointsData(0)->Fill((Int_t)clusters->GetEntriesFast()) ; 
534   
535   TIter next(clusters) ; 
536   AliTOFcluster * c ; 
537   while ( (c = dynamic_cast<AliTOFcluster *>(next())) ) {
538    
539       volumeID[0] = c->GetDetInd(0);
540       volumeID[1] = c->GetDetInd(1);
541       volumeID[2] = c->GetDetInd(2);
542       volumeID[3] = c->GetDetInd(4); //X and Z indeces inverted in RecPoints
543       volumeID[4] = c->GetDetInd(3); //X and Z indeces inverted in RecPoints
544       
545       for (Int_t i=0;i<5;i++){
546           volumeID2[i]=c->GetDetInd(i); //X and Z indeces inverted in RecPoints
547       }
548       
549       GetMapIndeces(volumeID,out);
550       Int_t chIndex=AliTOFGeometry::GetIndex(volumeID2);
551       Int_t iDDL=AliTOFRawStream::Geant2DDL(volumeID2);
552       Int_t iTRM=AliTOFRawStream::Geant2TRM(volumeID2);
553
554       if (c->GetTDCRAW() && c->GetTDC() && c->GetToT()){
555           if (volumeID2[0]>4 && volumeID2[0]<14){       //I side
556               if ((iDDL%4==0)|| (iDDL%4==1)){ //A side
557                   GetRecPointsData(1)->Fill( c->GetTDC()*tdc2ns) ;//in ns
558                   GetRecPointsData(5)->Fill( c->GetTDCRAW()*tdc2ns) ;//in ns
559                   GetRecPointsData(9)->Fill( c->GetToT()*tot2ns) ;//in ns
560                   
561               } else {
562                   if ((iDDL%4==2)|| (iDDL%4==3)){//C side
563                       GetRecPointsData(3)->Fill( c->GetTDC()*tdc2ns) ;//in ns
564                       GetRecPointsData(7)->Fill( c->GetTDCRAW()*tdc2ns) ;//in ns
565                       GetRecPointsData(11)->Fill( c->GetToT()*tot2ns) ;//in ns
566                   }
567               }
568           } else {
569               if (volumeID2[0]<5 || volumeID2[0]>13){       //O side
570                   if ((iDDL%4==0)|| (iDDL%4==1)){ //A side
571                       GetRecPointsData(2)->Fill( c->GetTDC()*tdc2ns) ;//in ns
572                       GetRecPointsData(6)->Fill( c->GetTDCRAW()*tdc2ns) ;//in ns
573                       GetRecPointsData(10)->Fill( c->GetToT()*tot2ns) ;//in ns
574                   } else {
575                       if ((iDDL%4==2)|| (iDDL%4==3)){//C side
576                           GetRecPointsData(4)->Fill( c->GetTDC()*tdc2ns) ;//in ns
577                           GetRecPointsData(8)->Fill( c->GetTDCRAW()*tdc2ns) ;//in ns
578                           GetRecPointsData(12)->Fill( c->GetToT()*tot2ns) ;//in ns
579                       }
580                   }
581               }
582           }
583           GetRecPointsData(14)->Fill(out[0],out[1]);
584           GetRecPointsData(15)->Fill(GetStripIndex(volumeID), c->GetTDC()*tdc2ns) ;
585           Int_t trm= iDDL*10+(iTRM-3);
586           if (iDDL>=0 && iDDL<36) {
587               GetRecPointsData(16)->Fill(trm,c->GetTDC()*tdc2ns);
588               GetRecPointsData(18)->Fill(trm,c->GetToT()*tot2ns);
589           }
590           if (iDDL>=36 && iDDL<72) {
591               GetRecPointsData(17)->Fill(trm,c->GetTDC()*tdc2ns);
592               GetRecPointsData(19)->Fill(trm,c->GetToT()*tot2ns);
593           }
594           GetRecPointsData(13)->Fill(chIndex) ;//in ns
595       }//hit selection
596   }//end while  
597  
598   EnableDqmShifterOpt(kFALSE);
599 }
600
601 //____________________________________________________________________________
602 void AliTOFQADataMakerRec::MakeESDs(AliESDEvent * esd)
603 {
604   //
605   // make QA data from ESDs
606   //  
607     const Double_t speedOfLight = TMath::C()*1E2*1E-12; // cm/ps
608     const Double_t pionMass = 0.13957018; //GeV/c^2
609
610     Int_t ntrk = esd->GetNumberOfTracks() ; 
611     Int_t ntof=0;
612     Int_t ntofpid=0;
613     Int_t nesdpid=0;
614     Int_t ntpc=0;
615     Int_t ntpctof=0;
616     Int_t ntofout=0;
617     
618     while (ntrk--) {
619         AliESDtrack *track=esd->GetTrack(ntrk);
620         Double_t tofTime=track->GetTOFsignal();//in ps
621         Double_t tofTimeRaw=track->GetTOFsignalRaw();//in ps
622         Double_t tofToT=track->GetTOFsignalToT(); //in ps
623
624         UInt_t status=track->GetStatus();
625         if (track->IsOn(AliESDtrack::kTPCrefit)) {
626             ntpc++;
627             Double_t y=track->Eta();
628             if ((status&AliESDtrack::kTOFout)!=0) {
629                 ntofout++;
630                 if (tofTime>0){
631                     ntof++;
632                     
633                     if (TMath::Abs(y)<0.9) {
634                         GetESDsData(6)->Fill(track->Pt());
635                         ntpctof++;
636                     }
637                     GetESDsData(1)->Fill(tofTime*1E-3);
638                     GetESDsData(2)->Fill(tofTimeRaw*1E-3); 
639                     GetESDsData(3)->Fill(tofToT*1E-3);
640                     //check how many tracks where ESD PID is ok 
641                     if ((status&AliESDtrack::kESDpid)!=0) nesdpid++; 
642                     if (((status&AliESDtrack::kESDpid)&AliESDtrack::kTOFpid)!=0) ntofpid++;
643                     
644                     Double_t length =track->GetIntegratedLength();
645                     Double_t mom2=(track->Pt()*track->Pt())+(track->Pz()*track->Pz());
646                     Double_t eTexp = TMath::Sqrt(1+(pionMass*pionMass/mom2))*length/speedOfLight; //in ps
647                     GetESDsData(9)->Fill(tofTime-eTexp);
648                     GetESDsData(10)->Fill(length);
649                 } //end check on matched tracks
650             }
651         }//end check on TPCrefit
652     }
653     
654     GetESDsData(0)->Fill(ntof) ;
655   
656     
657     if(ntof>0) {
658         Float_t ratio = (Int_t)ntofpid/(Int_t)ntof*100; //identified by TOF over matched
659         GetESDsData(4)->Fill(ratio) ;
660     }
661     
662     if(nesdpid>0) {
663         Float_t ratio = (Float_t)ntofpid/(Float_t)nesdpid *100; //identified by TOF over identified
664         GetESDsData(5)->Fill(ratio) ;
665     }
666     
667     if(ntpc>0){
668         Float_t ratio = (Float_t)ntof/(Float_t)ntpc*100.; //matching probability
669         GetESDsData(7)->Fill(ratio) ;
670     }
671     
672     if(ntofout>0) {
673         Float_t ratio = (Float_t)ntof/(Float_t)ntofout*100; //matched over propagated to TOF outer radius
674         GetESDsData(8)->Fill(ratio) ;
675     }
676     EnableDqmShifterOpt(kFALSE);
677 }
678
679 //____________________________________________________________________________ 
680 void AliTOFQADataMakerRec::StartOfDetectorCycle()
681 {
682   //
683   //Detector specific actions at start of cycle
684  
685   fCalibData = GetCalibData();
686
687 }
688
689 //____________________________________________________________________________ 
690 void AliTOFQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
691 {
692   //Detector specific actions at end of cycle
693   // do the QA checking
694     
695     for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
696         if ( !AliQAv1::Instance()->IsEventSpecieSet(specie) ) 
697             continue ; 
698         
699         AliInfo(Form("Processed %i physics raw events",fProcessedRawEventN));
700
701         if (fEnableDqmShifterOpt){
702             // Help make the raw qa histogram easier to interpret for the DQM shifter
703             // This is still to be optimized...
704
705             if (!GetRawsData(0) || !GetRawsData(5) || !GetRawsData(10) 
706                 || !GetRawsData(15) || !GetRawsData(16) || !GetRawsData(17)) {
707                 printf("No histogram for DQM found - Possible memory corruption ???. Please check\n") ; 
708                 continue;
709             }
710             
711             Double_t yTimeMin = GetRawsData(5)->GetMinimum();
712             Double_t yTimeMax = GetRawsData(5)->GetMaximum();
713             Double_t yTotMin = GetRawsData(10)->GetMinimum();
714             Double_t yTotMax = GetRawsData(10)->GetMaximum();
715             
716             TLine* lineExpTimeMin = new TLine(200., yTimeMin, 200., yTimeMax);
717             lineExpTimeMin->SetLineColor(kGreen);
718             lineExpTimeMin->SetLineWidth(2);
719             
720             TLine* lineExpTimeMax = new TLine(250., yTimeMin, 250., yTimeMax);
721             lineExpTimeMax->SetLineColor(kGreen);
722             lineExpTimeMax->SetLineWidth(2);
723             
724             TLine* lineExpTotMin = new TLine( 5., yTotMin, 5., yTotMax);
725             lineExpTotMin->SetLineColor(kGreen);
726             lineExpTotMin->SetLineWidth(2);
727             
728             TLine* lineExpTotMax = new TLine(20., yTotMin, 20., yTotMax);
729             lineExpTotMax->SetLineColor(kGreen);
730             lineExpTotMax->SetLineWidth(2);
731             
732             for(Int_t j=0;j<5;j++){
733                 //make up for hits count
734                 GetRawsData(j)->SetMarkerColor(kRed);
735                 GetRawsData(j)->SetMarkerStyle(7);
736                 //make up for time spectra
737                 ((TH1F*)GetRawsData(j+5))->GetListOfFunctions()->Add(lineExpTimeMin);
738                 ((TH1F*)GetRawsData(j+5))->GetListOfFunctions()->Add(lineExpTimeMax);
739                 GetRawsData(j+5)->SetLineColor(kBlue+1);
740                 GetRawsData(j+5)->SetMarkerColor(kBlue+1);
741                 
742                 //make up for tot spectra
743                 ((TH1F*)GetRawsData(j+10))->GetListOfFunctions()->Add(lineExpTotMin);
744                 ((TH1F*)GetRawsData(j+10))->GetListOfFunctions()->Add(lineExpTotMax);
745                 GetRawsData(j+10)->SetLineColor(kBlue+1);
746                 GetRawsData(j+10)->SetMarkerColor(kBlue+1);
747             }
748             
749             //make up for equipment hits count 
750             for(Int_t j=15;j<18;j++){
751               GetRawsData(j)->SetLineColor(kBlue+1);
752               GetRawsData(j)->SetLineWidth(1);
753
754               if (j==15) {
755                 GetRawsData(j)->SetMarkerStyle(8);
756                 GetRawsData(j)->SetMarkerSize(0.7);
757                 GetRawsData(j)->SetMarkerColor(kBlue+2);
758            
759               } else {
760                 Int_t ySMmax=GetRawsData(j)->GetMaximum();
761                 TLine* lineSMid[10];
762                 
763                 for (Int_t sm=0;sm<10;sm++){
764                   lineSMid[sm] = new TLine( 40*sm+360*(j%16), 0, 40*sm+360*(j%16), ySMmax);
765                   lineSMid[sm]->SetLineColor(kMagenta);
766                   lineSMid[sm]->SetLineWidth(2);
767                   GetRawsData(j)->GetListOfFunctions()->Add(lineSMid[sm]);
768                   GetRawsData(j)->SetFillColor(kBlue+1);
769                 }
770               }
771             }
772             GetRawsData(18)->SetLineColor(kBlue+1);
773             GetRawsData(18)->SetMarkerColor(kBlue+1);
774             
775             //make up for orphans histo
776             GetRawsData(19)->SetLineColor(kBlue+1);
777             GetRawsData(19)->SetMarkerColor(kBlue+1);
778             ((TH1F*)GetRawsData(19))->GetListOfFunctions()->Add(lineExpTimeMin);
779             ((TH1F*)GetRawsData(19))->GetListOfFunctions()->Add(lineExpTimeMax);
780         }
781     }
782     AliQAChecker::Instance()->Run(AliQAv1::kTOF, task, list) ;  
783 }
784 //____________________________________________________________________________
785 void AliTOFQADataMakerRec::GetMapIndeces(const Int_t* const in , Int_t* out)
786 {
787   //
788   //return appropriate indeces for the theta-phi map
789   //
790
791   Int_t npadX = AliTOFGeometry::NpadX();
792   Int_t npadZ = AliTOFGeometry::NpadZ();
793   Int_t nStripA = AliTOFGeometry::NStripA();
794   Int_t nStripB = AliTOFGeometry::NStripB();
795   Int_t nStripC = AliTOFGeometry::NStripC();
796
797   Int_t isector = in[0];
798   Int_t iplate = in[1];
799   Int_t istrip = in[2];
800   Int_t ipadX = in[3]; 
801   Int_t ipadZ = in[4]; 
802   
803   Int_t stripOffset = 0;
804   switch (iplate) {
805   case 0:
806     stripOffset = 0;
807       break;
808   case 1:
809     stripOffset = nStripC;
810     break;
811   case 2:
812     stripOffset = nStripC+nStripB;
813     break;
814   case 3:
815     stripOffset = nStripC+nStripB+nStripA;
816     break;
817   case 4:
818     stripOffset = nStripC+nStripB+nStripA+nStripB;
819     break;
820   default:
821     AliDebug(1,Form("Wrong plate number in TOF (%d) !",iplate));
822     break;
823   };
824   Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
825   Int_t phiindex=npadX*isector+ipadX+1;
826   out[0]=zindex;  
827   out[1]=phiindex;  
828   
829 }
830
831 //---------------------------------------------------------------
832 Int_t AliTOFQADataMakerRec::GetStripIndex(const Int_t * const in)
833 {
834     /* return tof strip index between 0 and 91 */
835
836   Int_t nStripA = AliTOFGeometry::NStripA();
837   Int_t nStripB = AliTOFGeometry::NStripB();
838   Int_t nStripC = AliTOFGeometry::NStripC();
839
840   // Int_t isector = in[0];
841   Int_t iplate = in[1];
842   Int_t istrip = in[2];
843   //Int_t ipadX = in[3]; 
844   //Int_t ipadZ = in[4]; 
845   
846   Int_t stripOffset = 0;
847   switch (iplate) {
848   case 0:
849     stripOffset = 0;
850       break;
851   case 1:
852     stripOffset = nStripC;
853     break;
854   case 2:
855     stripOffset = nStripC+nStripB;
856     break;
857   case 3:
858     stripOffset = nStripC+nStripB+nStripA;
859     break;
860   case 4:
861     stripOffset = nStripC+nStripB+nStripA+nStripB;
862     break;
863   default:
864       AliDebug(1,Form("Wrong plate number in TOF (%d) !",iplate));
865       stripOffset=-1;
866       break;
867   };
868   
869   if (stripOffset<0 || stripOffset>92) return -1;
870   else 
871       return (stripOffset+istrip);
872   
873 }
874 //---------------------------------------------------------------
875 Bool_t  AliTOFQADataMakerRec::CheckVolumeID(const Int_t * const volumeID)
876 {
877     //
878     //Checks volume ID validity
879     //
880    
881     for (Int_t j=0;j<5;j++){
882         if (volumeID[j]<0) {
883             AliDebug(1,Form("Invalid detector volume index for volumeID[%i]",j));
884             return kFALSE;
885         }
886     }
887     return kTRUE;
888     
889 }
890
891 //---------------------------------------------------------------
892 Bool_t  AliTOFQADataMakerRec::CheckEquipID(const Int_t * const equipmentID)
893 {
894     //
895     //Checks equipment ID validity
896     
897    for (Int_t j=0;j<5;j++){
898         if (equipmentID[j]<0) {
899           AliDebug(1,Form("Invalid equipment volume index for equipmentID[%i]",j));
900           return kFALSE;
901         }
902    }
903    return kTRUE;
904 }
905 //---------------------------------------------------------------
906 Bool_t  AliTOFQADataMakerRec::FilterLTMData(const Int_t * const equipmentID) const
907 {
908   /*It returns kTRUE if data come from LTM.
909     It thus filters trigger-related signals  */
910
911   Int_t ddl, trm, tdc;
912   //if (!CheckEquipID(equipmentID)) return kFALSE;
913   ddl = equipmentID[0];
914   trm = equipmentID[1];
915   tdc = equipmentID[3];
916   
917   if ((ddl%2==1) && (trm==3) && (tdc>11 && tdc<15))
918     return kTRUE;
919   else 
920     return kFALSE;
921  
922 }
923 //---------------------------------------------------------------
924 Bool_t  AliTOFQADataMakerRec::FilterSpare(const Int_t * const equipmentID) const
925 {
926   /*It returns kTRUE if data come from spare 
927     equipment ID. 
928     So far only check on TRM 3 crate left is implemented */
929
930   Int_t ddl, trm, tdc;
931   //if (!CheckEquipID(equipmentID)) return kFALSE;
932   ddl = equipmentID[0];
933   trm = equipmentID[1];
934   tdc = equipmentID[3];
935   
936   if ((ddl%2==1) && (trm==3) && (tdc>2 && tdc<12))
937     return kTRUE;
938   else 
939     return kFALSE;
940  
941 }