]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFQADataMakerRec.cxx
Adding useful print out in case of issues
[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 01/11/2011
26 - removed TLines as functions
27 - changed shifters plots for 2012 DQM
28
29 Modified by fbellini on 01/11/2011
30 - added histograms for LTM monitoring
31 - fix for coverity
32
33 Modified by fbellini on 17/10/2011
34 - fix for memory leak in constructor
35 - added methods to read histos ranges from config file in DQM
36 - added CTTM maps + relative methods to retrieve CTTM numbering
37 - removed obslete comments
38
39 Modified by fbellini & rshanoian on 06/07/2011
40 - changes for trigger classes implementation
41 - fRunNumber added as private member
42 - added time vs BCID plot
43
44 Modified by fbellini on 18/01/2011
45 - reduced histo binning to reduce size 
46 - added decoding errors plot
47 - added channel maps and options for DQM shifters
48 - new list of recPoints and ESDs plots
49 - removed hTrmChannels035 and hTrmChannels3671
50
51  Modified by bvonhall on 03/11/2010
52    - modified declaration of hTrmChannels035 and hTrmChannels3671 in EndOfDetectorCycle()
53      to prevent memory corruption
54    
55  Modified by adecaro on 18/10/2010
56    - fTOFRawStream object set as private member
57  
58 Modified by fbellini on 13/09/2010
59   - Set TLines as private members
60   - Set image flag for expert histos
61
62 Modified by fbellini on 14/06/2010
63   - Updated plots
64   - use LoadRawDataBuffersV2()
65
66   Modified by fbellini on 10/05/2010
67   - Fixed EndOfDetectorCycle() memory corruption bug
68
69   Modified by fbellini on 22/04/2010
70    - Added filter for physics events
71
72    Modified by fbellini on 16/04/2010
73    - Added EnableDqmShifterOpt() 
74    - Modified EndOfDetectorCycle() with options for DQM         
75    - Updated ESDs QA
76
77    Modified by fbellini on 30/03/2010
78    - Changed raws time histos range to 610ns
79    - Added FilterLTMData() and FilterSpare() methods
80    - Added check on enabled channels for raw data               
81    - Updated RecPoints QA
82
83    Modified by fbellini on 02/03/2010
84    - Fixed raw data decoding methods (use AliTOFRawStream::LoadRawDataBuffer())
85    - Added filter for noisy channels and read map from OCDB
86    - Added GetCalibData() method
87    - Added CheckVolumeID() and CheckEquipID() methods  
88    - Updated Raw QA
89 */
90 #include <iostream>
91 #include <fstream>
92
93 #include <TClonesArray.h>
94 #include <TH1F.h> 
95 #include <TH2F.h> 
96 #include <TLine.h>
97 #include <TPaveText.h>
98
99 #include "AliCDBManager.h"
100 #include "AliCDBEntry.h"
101 #include "AliESDEvent.h"
102 #include "AliESDtrack.h"
103 #include "AliQAChecker.h"
104 #include "AliRawReader.h"
105 #include "AliTOFRawStream.h"
106 #include "AliTOFcluster.h"
107 #include "AliTOFQADataMakerRec.h"
108 #include "AliTOFrawData.h"
109 #include "AliTOFGeometry.h"
110 #include "AliTOFChannelOnlineStatusArray.h"
111 #include "AliTOFDecoderSummaryData.h"
112 #include "AliTOFDecoderV2.h"
113
114 ClassImp(AliTOFQADataMakerRec)
115
116 Int_t AliTOFQADataMakerRec::fgNbinsMultiplicity=200; //number of bins in multiplicity plot
117 Int_t AliTOFQADataMakerRec::fgRangeMinMultiplicity=0;//min range in multiplicity plot
118 Int_t AliTOFQADataMakerRec::fgRangeMaxMultiplicity=200;//max range in multiplicity plot
119 Int_t AliTOFQADataMakerRec::fgNbinsTime=250;//number of bins in time plot
120 const Float_t AliTOFQADataMakerRec::fgkNbinsWidthTime=2.44;//width of bins in time plot
121 Float_t AliTOFQADataMakerRec::fgRangeMinTime=0.0;//range min in time plot
122 Float_t AliTOFQADataMakerRec::fgRangeMaxTime=620.0; //range max in time plot
123 Int_t AliTOFQADataMakerRec::fgCutNmaxFiredMacropad=50;//cut on number of max fired macropad
124 const Int_t AliTOFQADataMakerRec::fgkFiredMacropadLimit=50;//cut on number of max fired macropad
125
126
127 //____________________________________________________________________________ 
128   AliTOFQADataMakerRec::AliTOFQADataMakerRec() : 
129   AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kTOF), "TOF Quality Assurance Data Maker"),
130   fCalibData(0x0),
131   fEnableNoiseFiltering(kFALSE),
132   fEnableDqmShifterOpt(kFALSE),
133   fIsSOC(kFALSE),
134   fLineExpTimeMin(0x0),
135   fLineExpTimeMax(0x0),
136   fLineExpTotMin(0x0),
137   fLineExpTotMax(0x0),
138   fTOFRawStream(AliTOFRawStream()),
139   fDecoderSummary(0),
140   fRunNumber(-1),
141   fCalib(AliTOFcalib())
142 {
143   //
144   // ctor
145   //   
146   // fLineExpTimeMin = new TLine(200., 0., 200., 0.);
147   // fLineExpTimeMax = new TLine(250., 0., 250., 0.);
148   // fLineExpTotMin = new TLine(5., 0., 5., 0.);
149   // fLineExpTotMax = new TLine(20., 0., 20., 0.);
150   /*
151     for (Int_t sm=0;sm<17;sm++){
152     fLineSMid[sm] = new TLine( sm+1, 0., sm+1, 91.);
153   }
154
155   for (Int_t sm=0;sm<71;sm++){
156     fLineLTMid[sm] = new TLine( sm+1, 0., sm+1, 23.);
157   }
158
159   for (Int_t sm=0;sm<22;sm++){
160     fLineLTMbitId[sm] = new TLine( 0., sm+1, 72. ,sm+1);
161   }
162   */
163 }
164
165 //____________________________________________________________________________ 
166 AliTOFQADataMakerRec::AliTOFQADataMakerRec(const AliTOFQADataMakerRec& qadm) :
167   AliQADataMakerRec(),
168   fCalibData(qadm.fCalibData),
169   fEnableNoiseFiltering(qadm.fEnableNoiseFiltering),
170   fEnableDqmShifterOpt(qadm.fEnableDqmShifterOpt),
171   fIsSOC(qadm.fIsSOC),
172   fLineExpTimeMin(qadm.fLineExpTimeMin),
173   fLineExpTimeMax(qadm.fLineExpTimeMax),
174   fLineExpTotMin(qadm.fLineExpTotMin),
175   fLineExpTotMax(qadm.fLineExpTotMax),
176   fTOFRawStream(qadm.fTOFRawStream),
177   fDecoderSummary(qadm.fDecoderSummary),
178   fRunNumber(qadm.fRunNumber),
179   fCalib(qadm.fCalib)
180 {
181   //
182   //copy ctor 
183   //
184   SetName((const char*)qadm.GetName()) ; 
185   SetTitle((const char*)qadm.GetTitle()); 
186   /*
187   for (Int_t sm=0;sm<17;sm++){
188     fLineSMid[sm]=qadm.fLineSMid[sm];
189   }
190
191  for (Int_t sm=0;sm<71;sm++){
192     fLineLTMid[sm] = qadm.fLineLTMid[sm];
193   }
194
195   for (Int_t sm=0;sm<22;sm++){
196     fLineLTMbitId[sm] = qadm.fLineLTMbitId[sm];
197   }
198   */
199 }
200
201 //__________________________________________________________________
202 AliTOFQADataMakerRec& AliTOFQADataMakerRec::operator = (const AliTOFQADataMakerRec& qadm )
203 {
204   //
205   // assignment operator.
206   //
207   this->~AliTOFQADataMakerRec();
208   new(this) AliTOFQADataMakerRec(qadm);
209   return *this;
210 }
211  
212 //----------------------------------------------------------------------------
213 AliTOFQADataMakerRec::~AliTOFQADataMakerRec()
214 {
215   //destructor
216   fTOFRawStream.Clear();
217   fCalib.Clear();
218   if (fLineExpTimeMin)
219     delete fLineExpTimeMin;
220   if (fLineExpTimeMax)
221     delete fLineExpTimeMax;
222   if (fLineExpTotMin)
223     delete fLineExpTotMin;
224   if (fLineExpTotMax)
225     delete fLineExpTotMax;
226   /*
227     for (Int_t sm=0;sm<17;sm++){
228     if (fLineSMid[sm])
229       delete fLineSMid[sm];
230   }
231   for (Int_t sm=0;sm<71;sm++){
232     if (fLineLTMid[sm])
233       delete fLineLTMid[sm];
234   }
235 for (Int_t sm=0;sm<22;sm++){
236     if (fLineLTMbitId[sm])
237       delete fLineLTMbitId[sm];
238   }
239   */
240 }
241 //----------------------------------------------------------------------------
242 AliTOFChannelOnlineStatusArray* AliTOFQADataMakerRec::GetCalibData() 
243 {
244   //
245   // Retrive TOF calib objects from OCDB
246   //
247   AliCDBManager *man = AliCDBManager::Instance();
248   AliCDBEntry *cdbe=0;
249  
250   if (fRun<=0) fRunNumber=145288; //reference run from LHC11a
251   else fRunNumber=fRun;
252   
253   if (man->GetRun()!=fRunNumber){
254     fRunNumber=man->GetRun();
255     AliWarning(Form("Run number mismatch found: setting it to value from current AliCDBManager instance = %i", fRunNumber));
256   }
257   cdbe = man->Get("TOF/Calib/Status",fRunNumber);
258   
259   if(!cdbe){
260     // for DQM online
261     AliWarning("Load of calibration data from default (alien://) storage failed!");
262     printf("Calibration data will be loaded from local storage - ok if on DQM station!");
263     man->SetDefaultStorage("local:///local/cdb/");
264     cdbe = man->Get("TOF/Calib/Status",fRun);
265     
266     if(!cdbe){
267       AliWarning("Load of calibration data from local DQM machine storage failed!");
268       AliWarning("Calibration data will be loaded from local ($ALICE_ROOT) storage ");
269       man->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
270       cdbe = man->Get("TOF/Calib/Status",fRunNumber);
271     }
272   }
273   // Retrieval of data in directory TOF/Calib/Data:
274   AliTOFChannelOnlineStatusArray * array = 0;
275   if (cdbe) {
276     printf("======= OCDB object for TOF retrieved from run %i in %s\n",fRunNumber,cdbe->GetName());
277     array = (AliTOFChannelOnlineStatusArray *)cdbe->GetObject();
278   }
279   if (!array)  AliFatal("No calibration data from calibration database !");
280   
281   fCalib.Init(fRunNumber);
282   return array;
283 }
284
285 //____________________________________________________________________________ 
286 void AliTOFQADataMakerRec::InitRaws()
287 {
288   //
289   // create Raws histograms in Raws subdir
290   //
291   ReadHistogramRangeFromFile(gSystem->Getenv("TOFDQMHISTO_CONFIGFILE"));
292   
293   const Bool_t expert   = kTRUE ; 
294   const Bool_t saveCorr = kTRUE ; 
295   const Bool_t image    = kTRUE ; 
296
297   TH1I * h0 =  new TH1I("hTOFRaws","TOF raw hit multiplicity; TOF raw hits number; Events ",fgNbinsMultiplicity, fgRangeMinMultiplicity, fgRangeMaxMultiplicity);
298   TH1I * h1 =  new TH1I("hTOFRawsIA","TOF raw hit multiplicity - I/A side; TOF raw hits number;Events ",fgNbinsMultiplicity, fgRangeMinMultiplicity, fgRangeMaxMultiplicity);
299   TH1I * h2 =  new TH1I("hTOFRawsOA","TOF raw hit multiplicity - O/A side; TOF raw hits number;Events ",fgNbinsMultiplicity, fgRangeMinMultiplicity, fgRangeMaxMultiplicity);
300   TH1I * h3 =  new TH1I("hTOFRawsIC","TOF raw hit multiplicity - I/C side; TOF raw hits number;Events ",fgNbinsMultiplicity, fgRangeMinMultiplicity, fgRangeMaxMultiplicity);
301   TH1I * h4 =  new TH1I("hTOFRawsOC","TOF raw hit multiplicity - O/C side; TOF raw hits number;Events ",fgNbinsMultiplicity, fgRangeMinMultiplicity, fgRangeMaxMultiplicity);
302
303   TH1F * h5  = new TH1F("hTOFRawsTime", "TOF Raws - Hit time (ns);Measured Hit time [ns];Hits", fgNbinsTime,fgRangeMinTime,fgRangeMaxTime); 
304   TH1F * h6  = new TH1F("hTOFRawsTimeIA", "TOF Raws - Hit time (ns) - I/A side;Measured Hit time [ns];Hits", fgNbinsTime,fgRangeMinTime,fgRangeMaxTime); 
305   TH1F * h7  = new TH1F("hTOFRawsTimeOA", "TOF Raws - Hit time (ns) - O/A side;Measured Hit time [ns];Hits", fgNbinsTime,fgRangeMinTime,fgRangeMaxTime); 
306   TH1F * h8  = new TH1F("hTOFRawsTimeIC", "TOF Raws - Hit time (ns) - I/C side;Measured Hit time [ns];Hits", fgNbinsTime,fgRangeMinTime,fgRangeMaxTime); 
307   TH1F * h9  = new TH1F("hTOFRawsTimeOC", "TOF Raws - Hit time (ns) - O/C side;Measured Hit time [ns];Hits", fgNbinsTime,fgRangeMinTime,fgRangeMaxTime); 
308   
309   TH1F * h10  = new TH1F("hTOFRawsToT", "TOF Raws - Hit ToT (ns);Measured Hit ToT (ns);Hits", 100, 0., 48.8); 
310   
311   TH1F * h11  = new TH1F("hTOFRawsToTIA", "TOF Raws - Hit ToT (ns) - I/A side;Measured Hit ToT (ns);Hits", 100, 0., 48.8); 
312   TH1F * h12  = new TH1F("hTOFRawsToTOA", "TOF Raws - Hit ToT (ns) - O/A side;Measured Hit ToT (ns);Hits", 100, 0., 48.8); 
313   TH1F * h13  = new TH1F("hTOFRawsToTIC", "TOF Raws - Hit ToT (ns) - I/C side;Measured Hit ToT (ns);Hits", 100, 0., 48.8); 
314   TH1F * h14  = new TH1F("hTOFRawsToTOC", "TOF Raws - Hit ToT (ns) - O/C side;Measured Hit ToT (ns);Hits", 100, 0., 48.8); 
315   
316   TH1F * h15 = new TH1F("hTOFRawsLTMHits", "LTMs OR signals; Crate; Counts",  72, 0., 72.);
317   TH2F * h16  = new TH2F("hTOFrefMap", "TOF enabled channel reference map;sector;strip",  72, 0., 18., 91, 0., 91.);
318   TH2F * h17  = new TH2F("hTOFRawHitMap","TOF raw hit map;sector;strip", 72, 0., 18., 91, 0., 91.);
319   TH2I * h18 = new TH2I("hTOFDecodingErrors","Decoding error monitoring; DDL; Error ", 72, 0, 72, 13,1,14);
320   
321   h18->GetYaxis()->SetBinLabel(1,"DRM ");
322   h18->GetYaxis()->SetBinLabel(2,"LTM ");
323   h18->GetYaxis()->SetBinLabel(3,"TRM 3 ");
324   h18->GetYaxis()->SetBinLabel(4,"TRM 4");
325   h18->GetYaxis()->SetBinLabel(5,"TRM 5");
326   h18->GetYaxis()->SetBinLabel(6,"TRM 6");
327   h18->GetYaxis()->SetBinLabel(7,"TRM 7");
328   h18->GetYaxis()->SetBinLabel(8,"TRM 8");
329   h18->GetYaxis()->SetBinLabel(9,"TRM 9");
330   h18->GetYaxis()->SetBinLabel(10,"TRM 10");
331   h18->GetYaxis()->SetBinLabel(11,"TRM 11");
332   h18->GetYaxis()->SetBinLabel(12,"TRM 12");
333   h18->GetYaxis()->SetBinLabel(13,"recovered");
334   
335   TH1F * h19  = new TH1F("hTOFOrphansTime", "TOF Raws - Orphans time (ns);Measured Hit time [ns];Hits",fgNbinsTime,fgRangeMinTime,fgRangeMaxTime); 
336   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.,fgNbinsTime,fgRangeMinTime,fgRangeMaxTime);
337   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.,fgNbinsTime, fgRangeMinTime,fgRangeMaxTime);
338   TH2F * h22 = new TH2F("hTOFTimeVsStrip","TOF raw hit time vs. MRPC (along z axis); MRPC index along z axis; Raws TOF time (ns) ", 91,0.,91,fgNbinsTime,fgRangeMinTime,fgRangeMaxTime); 
339   TH2F * h23 = new TH2F("hTOFtimeVsBCID","TOF time vs BCID; BCID; time (ns) ", 3564, 0., 3564.,fgNbinsTime,fgRangeMinTime,fgRangeMaxTime);
340   TH2F * h24 = new TH2F("hTOFchannelEfficiencyMap","TOF channels (HWok && efficient && !noisy && !problematic);sector;strip",  72, 0., 18., 91, 0., 91.);
341   TH2F * h25 = new TH2F("hTOFhitsCTTM","Map of hit pads according to CTTM numbering;LTM index;bit index",  72, 0., 72., 23, 0., 23.);
342   TH2F * h26 = new TH2F("hTOFmacropadCTTM","Map of hit macropads according to CTTM numbering;LTM index; bit index",  72, 0., 72., 23, 0., 23.);
343   TH2F * h27 = new TH2F("hTOFmacropadDeltaPhiTime","#Deltat vs #Delta#Phi of hit macropads;#Delta#Phi (degrees);#DeltaBX",  18, 0., 180., 20, 0., 20.0);
344   TH2I *h28 = new TH2I("hBXVsCttmBit","BX ID in TOF matching window vs trg channel; trg channel; BX", 1728, 0, 1728, 24, 0, 24); 
345   TH2F *h29 = new TH2F("hTimeVsCttmBit","TOF raw time vs trg channel; trg channel; raw time (ns)", 1728, 0., 1728., fgNbinsTime, fgRangeMinTime, fgRangeMaxTime); 
346
347   h25->GetYaxis()->SetTickLength(-0.02);
348   h26->GetYaxis()->SetTickLength(-0.02);
349   h25->GetYaxis()->SetNdivisions(210);
350   h26->GetYaxis()->SetNdivisions(210);
351   h25->GetXaxis()->SetTickLength(-0.02);
352   h26->GetXaxis()->SetTickLength(-0.02);
353   h25->GetXaxis()->SetLabelOffset(0.015);
354   h26->GetXaxis()->SetLabelOffset(0.015);
355   h25->GetXaxis()->SetNdivisions(515);
356   h26->GetXaxis()->SetNdivisions(515);
357   
358   h0->Sumw2() ;
359   h1->Sumw2() ;
360   h2->Sumw2() ;
361   h3->Sumw2() ;
362   h4->Sumw2() ;
363   //  h5->Sumw2() ;
364   h6->Sumw2() ;
365   h7->Sumw2() ;
366   h8->Sumw2() ;
367   h9->Sumw2() ;
368   // h10->Sumw2() ;
369   h11->Sumw2() ;
370   h12->Sumw2() ;
371   h13->Sumw2() ;
372   h14->Sumw2() ;
373   h15->Sumw2() ;
374   h16->Sumw2() ;
375   h17->Sumw2() ;
376   h18->Sumw2() ;
377   h19->Sumw2() ;
378   h20->Sumw2() ;
379   h21->Sumw2() ;
380   h22->Sumw2() ;
381   h23->Sumw2() ;
382   h24->Sumw2() ;
383   h25->Sumw2() ;
384   h26->Sumw2() ;
385   h27->Sumw2() ;
386   h28->Sumw2() ;
387   h29->Sumw2() ;
388
389   //add lines for DQM shifter
390   fLineExpTimeMin = new TLine(200., 0., 200., 0.);
391   fLineExpTimeMax = new TLine(300., 0., 300., 0.);
392   fLineExpTotMin = new TLine(10., 0., 10., 0.);
393   fLineExpTotMax = new TLine(15., 0., 15., 0.);
394
395   fLineExpTimeMin->SetLineColor(kGreen);
396   fLineExpTimeMin->SetLineWidth(2);
397   
398   fLineExpTimeMax->SetLineColor(kGreen);
399   fLineExpTimeMax->SetLineWidth(2);
400   
401   fLineExpTotMin->SetLineColor(kGreen);
402   fLineExpTotMin->SetLineWidth(2);
403   
404   fLineExpTotMax->SetLineColor(kGreen);
405   fLineExpTotMax->SetLineWidth(2);
406   
407   /*
408     for (Int_t sm=0;sm<17;sm++){
409     fLineSMid[sm]->SetLineColor(kMagenta);
410     fLineSMid[sm]->SetLineWidth(2);
411   }
412   */
413   h5->GetListOfFunctions()->Add(fLineExpTimeMin);
414   h5->GetListOfFunctions()->Add(fLineExpTimeMax);
415   h10->GetListOfFunctions()->Add(fLineExpTotMin);
416   h10->GetListOfFunctions()->Add(fLineExpTotMax);
417   /*
418   for (Int_t sm=0;sm<17;sm++){
419     h16->GetListOfFunctions()->Add(fLineSMid[sm]);
420     h17->GetListOfFunctions()->Add(fLineSMid[sm]);
421   }
422   for (Int_t sm=0;sm<71;sm++){
423     fLineLTMid[sm]->SetLineColor(kBlack);
424     fLineLTMid[sm]->SetLineWidth(1);
425     h26->GetListOfFunctions()->Add(fLineLTMid[sm]);
426     h25->GetListOfFunctions()->Add(fLineLTMid[sm]);
427   }
428   for (Int_t sm=0;sm<22;sm++){
429     fLineLTMbitId[sm]->SetLineColor(kBlack);
430     fLineLTMbitId[sm]->SetLineWidth(1);
431     h26->GetListOfFunctions()->Add(fLineLTMbitId[sm]);
432     h25->GetListOfFunctions()->Add(fLineLTMbitId[sm]);
433   }
434   */
435   TPaveText *phosHoleBox=new TPaveText(13,38,16,53,"b");        
436   phosHoleBox->SetFillStyle(0);
437   phosHoleBox->SetFillColor(kWhite);
438   phosHoleBox->SetLineColor(kMagenta);
439   phosHoleBox->SetLineWidth(2);
440   phosHoleBox->AddText("PHOS"); 
441   h16->GetListOfFunctions()->Add(phosHoleBox);
442   h17->GetListOfFunctions()->Add(phosHoleBox);
443
444   // h0->SetDrawOption("logy");
445   // h5->SetDrawOption("logy");
446   // h10->SetDrawOption("logy");
447
448   Add2RawsList(h0,   0, !expert,  image, !saveCorr) ;
449   Add2RawsList(h1,   1,  expert, !image, !saveCorr) ;
450   Add2RawsList(h2,   2,  expert, !image, !saveCorr) ;
451   Add2RawsList(h3,   3,  expert, !image, !saveCorr) ;
452   Add2RawsList(h4,   4,  expert, !image, !saveCorr) ;
453   Add2RawsList(h5,   5, !expert,  image, !saveCorr) ;
454   Add2RawsList(h6,   6,  expert, !image, !saveCorr) ;
455   Add2RawsList(h7,   7,  expert, !image, !saveCorr) ;
456   Add2RawsList(h8,   8,  expert, !image, !saveCorr) ;
457   Add2RawsList(h9,   9,  expert, !image, !saveCorr) ;
458   Add2RawsList(h10, 10, !expert,  image, !saveCorr) ;
459   Add2RawsList(h11, 11,  expert, !image, !saveCorr) ;
460   Add2RawsList(h12, 12,  expert, !image, !saveCorr) ;
461   Add2RawsList(h13, 13,  expert, !image, !saveCorr) ;
462   Add2RawsList(h14, 14,  expert, !image, !saveCorr) ;
463   Add2RawsList(h15, 15,  expert, !image, !saveCorr) ;
464   Add2RawsList(h16, 16, !expert,  image, !saveCorr) ;
465   Add2RawsList(h17, 17, !expert,  image, !saveCorr) ;
466   Add2RawsList(h18, 18,  expert, !image, !saveCorr) ;
467   Add2RawsList(h19, 19,  expert, !image, !saveCorr) ;
468   Add2RawsList(h20, 20, !expert, image, !saveCorr) ;
469   Add2RawsList(h21, 21, !expert, image, !saveCorr) ;
470   Add2RawsList(h22, 22,  expert, !image, !saveCorr) ;
471   Add2RawsList(h23, 23,  expert, !image, !saveCorr) ;
472   Add2RawsList(h24, 24,  expert, !image, !saveCorr) ;
473   Add2RawsList(h25, 25,  expert, !image, !saveCorr) ;
474   Add2RawsList(h26, 26,  expert,  image, !saveCorr) ;
475   Add2RawsList(h27, 27,  expert, !image, !saveCorr) ;
476   Add2RawsList(h28, 28,  expert, !image, !saveCorr) ;
477   Add2RawsList(h29, 29,  expert, !image, !saveCorr) ;
478   
479 //
480   ClonePerTrigClass(AliQAv1::kRAWS); // this should be the last line
481 }
482
483 //____________________________________________________________________________ 
484 void AliTOFQADataMakerRec::InitRecPoints()
485 {
486   //
487   // create RecPoints histograms in RecPoints subdir
488   //
489   
490   const Bool_t expert   = kTRUE ; 
491   const Bool_t image    = kTRUE ; 
492
493   TH1I * h0 = new TH1I("hTOFRecPoints",    "TOF RecPoints multiplicity ; TOF RecPoints number;Events",200, 0, 200) ;
494
495   TH1F * h1 = new TH1F("hTOFRecPointsTimeIA", "RecPoints Time Spectrum in TOF (ns)- I/A side; Calibrated TOF time [ns];Events", 250, 0., 610.) ; 
496   TH1F * h2 = new TH1F("hTOFRecPointsTimeOA", "RecPoints Time Spectrum in TOF (ns)- O/A side; Calibrated TOF time [ns];Events", 250, 0., 610.) ;
497   TH1F * h3 = new TH1F("hTOFRecPointsTimeIC", "RecPoints Time Spectrum in TOF (ns)- I/C side; Calibrated TOF time [ns];Events", 250, 0., 610.) ;
498   TH1F * h4 = new TH1F("hTOFRecPointsTimeOC", "RecPoints Time Spectrum in TOF (ns)- O/C side; Calibrated TOF time [ns];Events", 250, 0., 610.) ;
499   
500   TH1F * h5  = new TH1F("hTOFRecPointsRawTimeIA", "RecPoints raw Time Spectrum in TOF (ns)-I/A side; Measured TOF time [ns];Hits", 250, 0., 610.) ; 
501   TH1F * h6  = new TH1F("hTOFRecPointsRawTimeOA", "RecPoints raw Time Spectrum in TOF (ns)-O/A side; Measured TOF time [ns];Hits", 250, 0., 610.) ; 
502   TH1F * h7  = new TH1F("hTOFRecPointsRawTimeIC", "RecPoints raw Time Spectrum in TOF (ns)-I/C side; Measured TOF time [ns];Hits", 250, 0., 610.) ; 
503   TH1F * h8  = new TH1F("hTOFRecPointsRawTimeOC", "RecPoints raw Time Spectrum in TOF (ns)-O/C side; Measured TOF time [ns];Hits", 250, 0., 610.) ; 
504  
505   TH1F * h9   = new TH1F("hTOFRecPointsToTIA", "RecPoints ToT Spectrum in TOF (ns)-I/A side; Measured TOT [ns];Hits", 100, 0., 48.8 ) ; 
506   TH1F * h10  = new TH1F("hTOFRecPointsToTOA", "RecPoints ToT Spectrum in TOF (ns)-O/A side; Measured TOT [ns];Hits", 100, 0., 48.8 ) ; 
507   TH1F * h11  = new TH1F("hTOFRecPointsToTIC", "RecPoints ToT Spectrum in TOF (ns)-I/C side; Measured TOT [ns];Hits", 100, 0., 48.8 ) ; 
508   TH1F * h12  = new TH1F("hTOFRecPointsToTOC", "RecPoints ToT Spectrum in TOF (ns)-O/C side; Measured TOT [ns];Hits", 100, 0., 48.8 ) ; 
509   
510   TH2F * h13 = new TH2F("hTOFRecPointsClusMap","RecPoints map; sector ;strip", 72, 0., 18., 91, 0., 91.) ; 
511   TH2F * h14 = new TH2F("hTOFRecPointsTimeVsStrip","RecPoints TOF time vs. strip (theta); Strip index; RecPoints TOF time (ns) ",91, 0., 91., 250, 0., 610.) ;
512   TH2F * h15 = 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) ;
513   TH2F * h16 = 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) ;
514
515   h0->Sumw2() ;
516   h1->Sumw2() ;
517   h2->Sumw2() ;
518   h3->Sumw2() ;
519   h4->Sumw2() ;
520   h5->Sumw2() ;
521   h6->Sumw2() ;
522   h7->Sumw2() ;
523   h8->Sumw2() ;
524   h9->Sumw2() ;
525   h10->Sumw2() ;
526   h11->Sumw2() ;
527   h12->Sumw2() ;
528   h13->Sumw2() ;
529   h14->Sumw2() ;
530   h15->Sumw2() ;
531   h16->Sumw2() ;
532    
533   Add2RecPointsList(h0, 0,   !expert,  image) ;
534   Add2RecPointsList(h1, 1,   !expert,  image) ;
535   Add2RecPointsList(h2, 2,   !expert,  image) ;
536   Add2RecPointsList(h3, 3,   !expert,  image) ;
537   Add2RecPointsList(h4, 4,   !expert,  image) ;
538   Add2RecPointsList(h5, 5,    expert,  !image) ;
539   Add2RecPointsList(h6, 6,    expert, !image) ;
540   Add2RecPointsList(h7, 7,    expert, !image) ;
541   Add2RecPointsList(h8, 8,    expert, !image) ;
542   Add2RecPointsList(h9, 9,   !expert, !image) ;
543   Add2RecPointsList(h10, 10, !expert, !image) ;
544   Add2RecPointsList(h11, 11, !expert, !image) ;
545   Add2RecPointsList(h12, 12, !expert, !image) ;
546   Add2RecPointsList(h13, 13,  expert, !image) ;
547   Add2RecPointsList(h14, 14,  expert, image) ;
548   Add2RecPointsList(h15, 15,  expert, !image) ;
549   Add2RecPointsList(h16, 16,  expert, !image) ;
550   //
551   ClonePerTrigClass(AliQAv1::kRECPOINTS); // this should be the last line
552 }
553
554 //____________________________________________________________________________ 
555 void AliTOFQADataMakerRec::InitESDs()
556 {
557   //
558   //create ESDs histograms in ESDs subdir
559   //
560
561   const Bool_t expert   = kTRUE ; 
562   const Bool_t image    = kTRUE ; 
563
564   TH1I * h0  = new TH1I("hTOFESDs", "Number of matched TOF tracks per event;Number of TOF matched ESD tracks;Counts", 200, 0, 200) ;  
565   TH1F * h1  = new TH1F("hTOFESDsTime", "Matched  ESDs tracks: TOF Time spectrum; Calibrated TOF time [ns];Counts", 250, 0., 610. ) ; 
566   TH1F * h2  = new TH1F("hTOFESDsRawTime", "Matched ESDs tracks: TOF raw Time spectrum;Measured TOF time [ns];Counts", 250, 0., 610.) ; 
567   TH1F * h3  = new TH1F("hTOFESDsToT", "Matched ESDs tracks: TOF ToT spectrum; ESDs ToT [ns];Counts",100, 0., 48.8) ; 
568   TH1F * h4  = new TH1F("hTOFESDskTOFOUT", "p_{T}  distribution of tracks with kTOFout; p_{T} (GeV/c);Counts", 50, 0.20, 5.00) ;  
569   TH1F * h5  = new TH1F("hTOFESDskTIME", "p_{T}  distribution of tracks with kTOFout && kTIME; p_{T} (GeV/c);Counts", 50, 0.20, 5.00) ;  
570   TH1F * h6  = new TH1F("hTOFESDsMatched", "p_{T} distribution of tracks with kTOFout && TOFtime>0; p_{T} (GeV/c);Counts", 50, 0.20, 5.00) ;  
571   TH1F * h7  = new TH1F("hTOFESDsMatchingProb", "TPC-TOF track-matching probability;TOF matching probability (%)  ;Counts",101, -1.0, 100) ;  
572   TH1F * h8  = new TH1F("hTOFESDsDiffTime", "ESDs t_{TOF}-t_{exp,pi} spectrum in TOF (ps); t_{TOF}-t_{exp,pi} [ps];Counts", 200, -2440., 2440.) ; 
573   TH1F * h9  = new TH1F("hTOFHitsLength", "Matched ESDs tracks: Length Spectrum; Track length [cm];Counts", 800, 0., 800) ; 
574
575   h0->Sumw2() ;
576   h1->Sumw2() ;
577   h2->Sumw2() ;
578   h3->Sumw2() ;
579   h4->Sumw2() ;
580   h5->Sumw2() ;
581   h6->Sumw2() ;
582   h7->Sumw2() ;
583   h8->Sumw2() ;
584   h9->Sumw2() ;
585
586   Add2ESDsList(h0, 0, !expert,  image) ;
587   Add2ESDsList(h1, 1, !expert,  image) ;
588   Add2ESDsList(h2, 2,  expert,  !image) ;
589   Add2ESDsList(h3, 3, !expert,  !image) ;
590   Add2ESDsList(h4, 4,  expert,  image) ;
591   Add2ESDsList(h5, 5,  expert,  image) ;
592   Add2ESDsList(h6, 6,  expert,  image) ; 
593   Add2ESDsList(h7, 7,  expert,  image) ; 
594   Add2ESDsList(h8, 8,  expert,  !image) ; 
595   Add2ESDsList(h9, 9, !expert,  !image) ;
596   //
597   ClonePerTrigClass(AliQAv1::kESDS); // this should be the last line 
598 }
599
600
601 //____________________________________________________________________________
602 void AliTOFQADataMakerRec::MakeRaws(AliRawReader* rawReader)
603 {
604   //
605   // makes data from Raws
606   //
607   // AliLog::SetClassDebugLevel("AliRawReader",0);
608   // AliLog::SetClassDebugLevel("AliTOFRawStream",0);
609   // AliLog::SetClassDebugLevel("AliTOFDecoderV2",0);
610   AliLog::SetGlobalLogLevel(AliLog::kError);
611   if (rawReader->GetType()==7) {
612    
613     Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;//in ns
614     Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
615     Int_t ntof[5]; /* 0=tot, 1=IA, 2=OA, 3=IC, 4=OC*/
616     for (Int_t j=0;j<5;j++){ ntof[j]=0;}
617     Int_t equipmentID[5]; //(ddl, trm, chain,tdc,channel)
618     Int_t volumeID[5];   //(sector,plate,strip,padX,padZ)
619     Int_t volumeID2[5];   //(sector,plate,strip,padZ,padX) to use AliTOFGeometry::GetIndex()
620     Int_t chIndex=-1;
621     Int_t indexCTTM[2]={-1,-1};     
622     Int_t indexGeo2CTTM[2]={-1,-1};
623     Float_t macropadPhiTimeUPC[fgkFiredMacropadLimit][2];
624     for (Int_t ii=0;ii<2;ii++){
625       for (Int_t jj=0;jj<fgkFiredMacropadLimit;jj++){   
626         macropadPhiTimeUPC[jj][ii]=-999.0; 
627       }
628     }
629
630     TClonesArray * clonesRawData;
631     fTOFRawStream.SetRawReader(rawReader);
632     Int_t BCID=rawReader->GetBCID();
633     
634     Int_t nFiredMacropad=0,
635       iFiredMacropad=-1;
636     nFiredMacropad=GetNumberOfFiredMacropad(rawReader);
637     
638     //uncomment if needed to apply DeltaBC correction
639     //fTOFRawStream.ApplyBCCorrections(kTRUE);
640     
641     if (fDecoderSummary){
642       fDecoderSummary->Reset();
643     }
644     for (Int_t iDDL = 0; iDDL < AliTOFGeometry::NDDL()*AliTOFGeometry::NSectors(); iDDL++){
645       rawReader->Reset();
646       fTOFRawStream.LoadRawDataBuffersV2(iDDL);
647       
648       //* get decoding error counters
649       fDecoderSummary = ( (AliTOFDecoderV2*) fTOFRawStream.GetDecoderV2() )->GetDecoderSummaryData();
650       if ( (fDecoderSummary) && (fDecoderSummary ->GetErrorDetected()) ) {
651         Int_t errorSlotID=(Int_t) fDecoderSummary->GetErrorSlotID();
652         FillRawsData(18,iDDL,errorSlotID);
653         if (fDecoderSummary -> GetRecoverError() )              
654           FillRawsData(18,iDDL,13);
655       }     
656       
657       clonesRawData = (TClonesArray*)fTOFRawStream.GetRawData();
658       for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
659         AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
660         Float_t tofRawTime=tofRawDatum->GetTOF()*tdc2ns;
661         
662         if (tofRawDatum->GetTOF()){
663           equipmentID[0]=iDDL;
664           equipmentID[1]=tofRawDatum->GetTRM(); 
665           equipmentID[2]=tofRawDatum->GetTRMchain();
666           equipmentID[3]=tofRawDatum->GetTDC();
667           equipmentID[4]=tofRawDatum->GetTDCchannel();
668           
669           if (CheckEquipID(equipmentID)){
670             fTOFRawStream.EquipmentId2VolumeId(iDDL, 
671                                                tofRawDatum->GetTRM(), 
672                                                tofRawDatum->GetTRMchain(),
673                                                tofRawDatum->GetTDC(), 
674                                                tofRawDatum->GetTDCchannel(), 
675                                                volumeID);
676             //LTM data
677             if (FilterLTMData(equipmentID)) { //counts LTM hits
678               if (equipmentID[2]==1)  //crate left, A-side or C-side
679                 FillRawsData(15,equipmentID[0]);
680               else 
681                 FillRawsData(15,equipmentID[0]-1); 
682               
683               //retrieve CTTM index (Ltm, bit)
684               GetCTTMIndex(equipmentID, indexCTTM);
685               
686               //get BX index within TOF-matching window
687               Int_t indexBC=-1;
688               indexBC= TMath::Nint(tofRawTime/24.4);
689
690               Int_t indexCttmChannel=indexCTTM[0]*24+indexCTTM[1];
691               FillRawsData(28,indexCttmChannel,indexBC);
692               FillRawsData(29,indexCttmChannel,tofRawTime);
693               
694               //fired macropad map (from LTM hits) - only for low multi evts (UPC)
695               if ((nFiredMacropad<=fgCutNmaxFiredMacropad)){    
696                 iFiredMacropad++;
697                 //AliInfo(Form("Event found with %i fired macropads in BCID = %i!",nFiredMacropad,BCID));
698                 FillRawsData(26,indexCTTM[0],indexCTTM[1]);
699                 Float_t halfSMphi=-999.0;
700                 if (indexCTTM[0]<36)
701                   halfSMphi=indexCTTM[0]*10.+5.;
702                 else  halfSMphi=(indexCTTM[0]-36)*10.+5.;
703                 macropadPhiTimeUPC[iFiredMacropad][0]=halfSMphi;
704                 macropadPhiTimeUPC[iFiredMacropad][1]=indexBC;
705               }
706             }
707             
708             //TRM data
709             if (CheckVolumeID(volumeID)){  
710               volumeID2[0]=volumeID[0];
711               volumeID2[1]=volumeID[1];
712               volumeID2[2]=volumeID[2];
713               volumeID2[3]=volumeID[4];
714               volumeID2[4]=volumeID[3];
715               chIndex=AliTOFGeometry::GetIndex(volumeID2);
716               
717               //fill hit map according to CTTM numbering
718               GetGeo2CTTMIndex(volumeID2, indexGeo2CTTM);
719               if ((nFiredMacropad<=fgCutNmaxFiredMacropad)){    
720                 FillRawsData(25,indexGeo2CTTM[0],indexGeo2CTTM[1]);
721               }
722               //hits selection
723               if (tofRawDatum->GetTOT()){           
724                 if (!(fCalibData->GetNoiseStatus(chIndex)==AliTOFChannelOnlineStatusArray::kTOFNoiseBad)
725                     && (fCalibData->GetHWStatus(chIndex) == AliTOFChannelOnlineStatusArray::kTOFHWOk)) {//noise and enabled filter
726                   ntof[0]++; //counter for tof hits
727                   
728                   //fill global spectra for DQM plots
729                   FillRawsData(5, tofRawTime) ;//in ns
730                   FillRawsData(10, tofRawDatum->GetTOT()*tot2ns) ;//in ns
731                   FillRawsData(23, BCID, tofRawTime) ;//in ns
732                   
733                   //fill side-related spectra for experts plots
734                   Int_t ddlACside=iDDL/36; // 0 or 1
735                   Int_t ddlPerSm=iDDL%4;
736                   
737                   if (volumeID2[0]>4 && volumeID2[0]<14){       //O side
738                     if (ddlPerSm<2){ //A side
739                       ntof[1]++;
740                       FillRawsData(6, tofRawTime) ;
741                       FillRawsData(11, tofRawDatum->GetTOT()*tot2ns) ;
742                     } else {  //C side
743                       ntof[3]++;
744                       FillRawsData(8, tofRawTime) ;
745                       FillRawsData(13, tofRawDatum->GetTOT()*tot2ns) ;
746                     }
747                   } else {                                    
748                     if (volumeID2[0]<5 || volumeID2[0]>13){   //I side
749                       if (ddlPerSm<2){ //A side
750                         ntof[2]++;
751                         FillRawsData(7, tofRawTime) ;
752                         FillRawsData(12, tofRawDatum->GetTOT()*tot2ns) ;
753                       } else {//C side
754                         ntof[4]++;
755                         FillRawsData(9, tofRawTime) ;
756                         FillRawsData(14, tofRawDatum->GetTOT()*tot2ns) ;
757                       }
758                     }   
759                   }
760                   
761                   //compute TRM offset
762                   Int_t trm= iDDL*10+(equipmentID[1]-3);
763                   FillRawsData(20+ddlACside,trm,tofRawTime);
764                   FillRawsData(22,GetStripIndex(volumeID),tofRawTime) ;
765                   Short_t fea = volumeID2[4]/12;
766                   Float_t hitmapx = volumeID2[0] + ((Double_t)(3 - fea) + 0.5) *0.25;
767                   FillRawsData(17,hitmapx,GetStripIndex(volumeID2));
768                   if (fCalib.IsChannelEnabled(chIndex,kTRUE,kTRUE))//checks also if efficient and if problematic
769                     FillRawsData(24,hitmapx,GetStripIndex(volumeID2));
770                 }//noise filter
771               }//end hit selection
772               else { //orphans
773                 if (!(fCalibData->GetNoiseStatus(chIndex) == AliTOFChannelOnlineStatusArray::kTOFNoiseBad)
774                     && (fCalibData->GetHWStatus(chIndex) == AliTOFChannelOnlineStatusArray::kTOFHWOk))
775                   FillRawsData(19, tofRawTime) ;//in ns
776               }//end orphans
777             }//end volumeID check
778           }//end equipID check
779         }//end tof check
780       }//loop on raw data
781       clonesRawData->Clear();
782     } // DDL Loop
783     
784     for (Int_t j=0;j<5;j++) FillRawsData(j,ntof[j]);
785     fTOFRawStream.Clear();
786   
787     if ((nFiredMacropad<=fgCutNmaxFiredMacropad)){
788       Float_t deltaPhiMacropad=-999.;
789       Float_t deltaTimeMacropad=-999.;
790       for (Int_t j=0;j<fgCutNmaxFiredMacropad+1; j++){
791         for (Int_t k=j+1;k<fgCutNmaxFiredMacropad+1; k++){
792           if ((macropadPhiTimeUPC[j][0]>0.0)&&(macropadPhiTimeUPC[k][0]>0.0)){
793             deltaPhiMacropad=TMath::Abs(macropadPhiTimeUPC[j][0]-macropadPhiTimeUPC[k][0]);
794             deltaTimeMacropad=TMath::Abs(macropadPhiTimeUPC[j][1]-macropadPhiTimeUPC[k][1]);
795             if (deltaPhiMacropad<=180.)
796               FillRawsData(27, deltaPhiMacropad,deltaTimeMacropad);
797             else
798               FillRawsData(27, TMath::Abs(360.0-deltaPhiMacropad),deltaTimeMacropad);
799           }
800         } 
801       }    
802     }//end cut on number of fired macropad
803   } else {
804     AliDebug(2,Form("Event of type %d found. Skipping non-physics event for QA.\n", rawReader->GetType())); 
805   }
806   
807   //fill reference map for DQM shifter only once in a detector cycle 
808   if (fIsSOC) {
809     Int_t geoId[5]={-1,-1,-1,-1,-1};// pgeoId=(sm, mod, strip, padZ, padX)
810     Int_t detId[5]={-1,-1,-1,-1,-1};//detID=(ddl,trm,tdc, chain,channel)
811     Int_t trmIndex=-1;
812     for (Int_t ch = 0; ch <  fCalibData->GetSize(); ch++) {
813       AliTOFGeometry::GetVolumeIndices(ch,geoId);
814       AliTOFRawStream::Geant2EquipmentId(geoId,detId); 
815       if ((detId[1]<0)||(detId[0]<0)) continue;
816       trmIndex=(detId[1]-3)+detId[0]*10;
817
818       if ( (!(fCalibData->GetNoiseStatus(ch)==AliTOFChannelOnlineStatusArray::kTOFNoiseBad))
819            && (fCalibData->GetHWStatus(ch) == AliTOFChannelOnlineStatusArray::kTOFHWOk) ){
820         //fill reference map with info from OCDB
821         Short_t fea = geoId[4]/12;
822         Float_t hitmapx = geoId[0] + ((Double_t)(3 - fea) + 0.5)*0.25;
823         FillRawsData(16,hitmapx, GetStripIndex(geoId));
824       }
825     }
826     fIsSOC=kFALSE;
827   }
828   
829   //enable options for DQM shifter
830   EnableDqmShifterOpt(kTRUE);
831   //
832   IncEvCountCycleRaws();
833   IncEvCountTotalRaws();
834   //
835 }
836
837 //____________________________________________________________________________
838 void AliTOFQADataMakerRec::MakeRecPoints(TTree * clustersTree)
839 {
840     //
841   // Make data from Clusters
842   //
843  
844   Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
845   Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
846   
847   Int_t volumeID2[5];//(sm, plate,strip, padZ,padX)
848   //Int_t volumeID[5];//(sm, plate,strip, padX,padZ) 
849     
850   TBranch *branch=clustersTree->GetBranch("TOF");
851   if (!branch) { 
852     AliError("can't get the branch with the TOF clusters !");
853     return;
854   }
855
856   static TClonesArray dummy("AliTOFcluster",10000);
857   dummy.Clear();
858   TClonesArray *clusters=&dummy;
859   branch->SetAddress(&clusters);
860   
861   // Import the tree
862   clustersTree->GetEvent(0);  
863  
864   FillRecPointsData(0,(Int_t)clusters->GetEntriesFast()) ; 
865   
866   TIter next(clusters) ; 
867   AliTOFcluster * c ; 
868   while ( (c = dynamic_cast<AliTOFcluster *>(next())) ) {
869    
870       // volumeID2[0] = c->GetDetInd(0);
871       // volumeID2[1] = c->GetDetInd(1);
872       // volumeID2[2] = c->GetDetInd(2);
873       // volumeID2[3] = c->GetDetInd(4); //padX
874       // volumeID2[4] = c->GetDetInd(3); //padZ 
875       
876       for (Int_t i=0;i<5;i++){
877         volumeID2[i]=c->GetDetInd(i); //X and Z indeces inverted in RecPoints
878       }
879       //Int_t chIndex=AliTOFGeometry::GetIndex(volumeID2);
880       Int_t iDDL=AliTOFRawStream::Geant2DDL(volumeID2);
881       Int_t iTRM=AliTOFRawStream::Geant2TRM(volumeID2);
882       Short_t fea = volumeID2[4]/12;
883       Float_t hitmapx = volumeID2[0] + ((Double_t)(3 - fea) + 0.5) *0.25;
884       Int_t ddlACside=iDDL/36; // 0 or 1
885       Int_t ddlPerSm=iDDL%4;
886       
887       if ((c->GetTDCRAW()) && (c->GetTDC()) && (c->GetToT())){
888         if (volumeID2[0]>4 && volumeID2[0]<14){       //I side
889           if (ddlPerSm<2){ //A side
890             FillRecPointsData(1, c->GetTDC()*tdc2ns) ;//in ns
891             FillRecPointsData(5, c->GetTDCRAW()*tdc2ns) ;//in ns
892             FillRecPointsData(9, c->GetToT()*tot2ns) ;//in ns
893           } else {//C side
894             FillRecPointsData(3, c->GetTDC()*tdc2ns) ;//in ns
895             FillRecPointsData(7, c->GetTDCRAW()*tdc2ns) ;//in ns
896             FillRecPointsData(11, c->GetToT()*tot2ns) ;//in ns
897           }
898         } else {
899           if (volumeID2[0]<5 || volumeID2[0]>13){       //O side
900             if (ddlPerSm<2){ //A side
901               FillRecPointsData(2, c->GetTDC()*tdc2ns) ;//in ns
902               FillRecPointsData(6, c->GetTDCRAW()*tdc2ns) ;//in ns
903               FillRecPointsData(10, c->GetToT()*tot2ns) ;//in ns
904             } else { //C side
905               FillRecPointsData(4, c->GetTDC()*tdc2ns) ;//in ns
906               FillRecPointsData(8, c->GetTDCRAW()*tdc2ns) ;//in ns
907               FillRecPointsData(12, c->GetToT()*tot2ns) ;//in ns
908             }
909           }
910         }
911         FillRecPointsData(13,hitmapx,GetStripIndex(volumeID2));
912         FillRecPointsData(14,GetStripIndex(volumeID2), c->GetTDC()*tdc2ns) ;
913         Int_t trm= iDDL*10+(iTRM-3);
914         if (ddlACside==0) { //A side
915           FillRecPointsData(15,trm,c->GetTDC()*tdc2ns);
916         } else {//C side
917           FillRecPointsData(16,trm,c->GetTDC()*tdc2ns);
918         }
919       }//hit selection
920   }//end while   
921   EnableDqmShifterOpt(kFALSE);
922   //
923   IncEvCountCycleRecPoints();
924   IncEvCountTotalRecPoints();
925   //
926 }
927
928 //____________________________________________________________________________
929 void AliTOFQADataMakerRec::MakeESDs(AliESDEvent * esd)
930 {
931   //
932   // make QA data from ESDs
933   //  
934   const Double_t speedOfLight = TMath::C()*1E2*1E-12; // cm/ps
935   const Double_t pionMass = 0.13957018; //GeV/c^2
936
937   Int_t ntrk = esd->GetNumberOfTracks() ; 
938   Int_t ntpc=0;
939   Int_t ntofout=0;
940   
941     while (ntrk--) {
942       AliESDtrack *track=esd->GetTrack(ntrk);
943       Double_t tofTime=track->GetTOFsignal();//in ps
944       Double_t tofTimeRaw=track->GetTOFsignalRaw();//in ps
945       Double_t tofToT=track->GetTOFsignalToT(); //in ps
946       
947       UInt_t status=track->GetStatus();
948       if (track->IsOn(AliESDtrack::kTPCrefit)) {
949         ntpc++;
950         Double_t y=track->Eta();
951         if (TMath::Abs(y)<0.9) { //select TOF acceptance
952           if ((status&AliESDtrack::kTOFout)!=0)  { //define matching
953             ntofout++;
954             FillESDsData(1,tofTime*1E-3);
955             FillESDsData(2,tofTimeRaw*1E-3); 
956             FillESDsData(3,tofToT*1E-3);
957             FillESDsData(4,track->Pt());
958             
959             Double_t length =track->GetIntegratedLength();
960             Double_t mom2=(track->Pt()*track->Pt())+(track->Pz()*track->Pz());
961             Double_t piTexp = TMath::Sqrt(1+(pionMass*pionMass/mom2))*length/speedOfLight; //in ps
962             FillESDsData(8,tofTime-piTexp);
963             FillESDsData(9,length);
964             
965             if ((status&AliESDtrack::kTIME)!=0) 
966               FillESDsData(5,track->Pt());
967             
968             if (tofTime>0)
969               FillESDsData(6,track->Pt());
970           } //end check on matched tracks
971         } 
972       }//end check on TPCrefit
973     }
974     
975     FillESDsData(0,ntofout) ;
976     if(ntpc>0){
977       Float_t ratio = (Float_t)ntofout/(Float_t)ntpc*100.; //matching probability
978       FillESDsData(7,ratio) ;
979     }
980     
981     if(ntofout>0) {
982         Float_t ratio = (Float_t)ntofout/(Float_t)ntpc*100; //matched over propagated to TOF outer radius
983         FillESDsData(8,ratio) ;
984     }
985     EnableDqmShifterOpt(kFALSE);
986     //
987     IncEvCountCycleESDs();
988     IncEvCountTotalESDs();
989     //
990 }
991
992 //____________________________________________________________________________ 
993 void AliTOFQADataMakerRec::StartOfDetectorCycle()
994 {
995   //
996   //Detector specific actions at start of cycle
997   fCalibData = GetCalibData();
998   fIsSOC=kTRUE;
999   return;
1000 }  
1001
1002 //____________________________________________________________________________ 
1003 void AliTOFQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
1004 {
1005   //Detector specific actions at end of cycle
1006   // do the QA checking
1007   ResetEventTrigClasses();
1008   //
1009   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
1010     if ( !AliQAv1::Instance()->IsEventSpecieSet(specie) ) continue ;     
1011     SetEventSpecie(AliRecoParam::ConvertIndex(specie));  
1012
1013     for (int itc=-1;itc<GetNTrigClasses();itc++) { // RS: loop over eventual clones per trigger class
1014
1015       if (fEnableDqmShifterOpt) {
1016         // RS: fetch the histograms for given trigger class
1017         TObjArray& arrRW = *GetRawsDataOfTrigClass(itc);
1018         
1019         // Help make the raw qa histogram easier to interpret for the DQM shifter
1020         if (!arrRW[ 0] || !arrRW[ 5] || !arrRW[10] || !arrRW[15] || !arrRW[16] || !arrRW[17]) continue;
1021         
1022         printf("=========>Processed %i physics raw of specie %s with TrigGlass %d\n",
1023                GetEvCountCycleRaws(itc),AliRecoParam::GetEventSpecieName(specie), itc);
1024         
1025         //Double_t monitorPeriodLength=fProcessedRawEventN*600*1E-9;//in s
1026       
1027         if (fCalibData){
1028           //set minima and maxima to allow log scale
1029           Double_t yTimeMax = ((TH1*)arrRW[5])->GetMaximum()*1.05;
1030           Double_t yTotMax = ((TH1*)arrRW[10])->GetMaximum()*1.05;
1031           fLineExpTimeMin->SetY2(yTimeMax);
1032           fLineExpTimeMax->SetY2(yTimeMax);
1033           fLineExpTotMin->SetY2(yTotMax);
1034           fLineExpTotMax->SetY2(yTotMax);
1035           //
1036           for (Int_t j=0;j<18;j++){
1037             if ((j==0)||(j==5)||(j==10)||(j==15)||(j==16)||(j==17)) {
1038               TH1* htmp = (TH1*)arrRW[j];
1039               htmp->GetXaxis()->SetLabelOffset(0.005);
1040               htmp->GetXaxis()->SetLabelSize(0.05);
1041               htmp->GetXaxis()->SetTitleOffset(0.8);
1042               htmp->GetXaxis()->SetTitleSize(0.05);
1043               htmp->GetYaxis()->SetLabelOffset(0.005);
1044               htmp->GetYaxis()->SetLabelSize(0.06);
1045               htmp->GetYaxis()->SetTitleOffset(0.8);
1046               htmp->GetYaxis()->SetTitleSize(0.06);       
1047             }
1048           }
1049           //make up for all histos 
1050           for(Int_t j=0;j<5;j++) {
1051             TH1* htmp = (TH1*)arrRW[j];  
1052             if (!htmp) continue;
1053             htmp->SetMarkerColor(kBlue);
1054             htmp->SetMarkerStyle(8);
1055             htmp->SetMarkerSize(0.7);
1056           }
1057           for(Int_t j=5;j<15;j++) {
1058             TH1* htmp = (TH1*)arrRW[j];
1059             if (!htmp) continue;
1060             htmp->SetLineColor(kBlue);
1061             htmp->SetLineWidth(1);
1062             htmp->SetMarkerColor(kBlue);
1063           }
1064           
1065           TH1* htmp =  (TH1*)arrRW[15];  
1066           htmp->SetLineColor(kBlue);
1067           htmp->SetLineWidth(1);
1068           htmp->SetMarkerStyle(8);
1069           htmp->SetMarkerSize(0.7);
1070           htmp->SetMarkerColor(kBlue);//Option("bar");
1071           //
1072           TString title25 = Form("Map of hit pads according to CTTM numbering (Max Fired Macropad = %i)",fgCutNmaxFiredMacropad);
1073           TString title26 = Form("Map of hit macropads according to CTTM numbering (Max Fired Macropad = %i)",fgCutNmaxFiredMacropad);
1074           TString title27 = Form("#Deltat vs #Delta#Phi of hit macropads (Max Fired Macropad = %i)",fgCutNmaxFiredMacropad);
1075           
1076           if ( (htmp=(TH1*)arrRW[16]) ) htmp->SetOption("colz");
1077           if ( (htmp=(TH1*)arrRW[17]) ) htmp->SetOption("colz");
1078           if ( (htmp=(TH1*)arrRW[18]) ) htmp->SetOption("colz"); 
1079           if ( (htmp=(TH1*)arrRW[22]) ) htmp->SetOption("colz"); 
1080           if ( (htmp=(TH1*)arrRW[23]) ) htmp->SetOption("colz"); 
1081           if ( (htmp=(TH1*)arrRW[24]) ) htmp->SetOption("colz"); 
1082           if ( (htmp=(TH1*)arrRW[28]) ) htmp->SetOption("colz"); 
1083           if ( (htmp=(TH1*)arrRW[29]) ) htmp->SetOption("colz"); 
1084
1085           if ( (htmp=(TH1*)arrRW[25]) ) {
1086             htmp->SetOption("colz"); 
1087             htmp->SetTitle(title25.Data());
1088           }
1089           if ( (htmp=(TH1*)arrRW[26]) ) {
1090             htmp->SetOption("colz"); 
1091             htmp->SetTitle(title26.Data());
1092           }
1093           if ( (htmp=(TH1*)arrRW[27]) ){
1094             htmp->SetOption("colz"); 
1095             htmp->SetTitle(title27.Data());
1096           }
1097         }
1098       }//END ENABLE DQM SHIFTER OPT
1099     } // RS: loop over trigger classes
1100   } //end for
1101   //
1102   AliQAChecker::Instance()->Run(AliQAv1::kTOF, task, list) ;  
1103   //
1104 }
1105 //____________________________________________________________________________
1106 void AliTOFQADataMakerRec::GetMapIndeces(const Int_t* const in , Int_t* out)
1107 {
1108   //
1109   //return appropriate indeces for the theta-phi map
1110   //
1111
1112   Int_t npadX = AliTOFGeometry::NpadX();
1113   Int_t npadZ = AliTOFGeometry::NpadZ();
1114   Int_t nStripA = AliTOFGeometry::NStripA();
1115   Int_t nStripB = AliTOFGeometry::NStripB();
1116   Int_t nStripC = AliTOFGeometry::NStripC();
1117
1118   Int_t isector = in[0];
1119   Int_t iplate = in[1];
1120   Int_t istrip = in[2];
1121   Int_t ipadX = in[3]; 
1122   Int_t ipadZ = in[4]; 
1123   
1124   Int_t stripOffset = 0;
1125   switch (iplate) {
1126   case 0:
1127     stripOffset = 0;
1128       break;
1129   case 1:
1130     stripOffset = nStripC;
1131     break;
1132   case 2:
1133     stripOffset = nStripC+nStripB;
1134     break;
1135   case 3:
1136     stripOffset = nStripC+nStripB+nStripA;
1137     break;
1138   case 4:
1139     stripOffset = nStripC+nStripB+nStripA+nStripB;
1140     break;
1141   default:
1142     //    AliDebug(2,Form("Wrong plate number in TOF (%d) !",iplate));
1143     break;
1144   };
1145   Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
1146   Int_t phiindex=npadX*isector+ipadX+1;
1147   out[0]=zindex;  
1148   out[1]=phiindex;  
1149   
1150 }
1151
1152 //---------------------------------------------------------------
1153 Int_t AliTOFQADataMakerRec::GetStripIndex(const Int_t * const in)
1154 {
1155     /* return tof strip index between 0 and 91 */
1156
1157   Int_t nStripA = AliTOFGeometry::NStripA();
1158   Int_t nStripB = AliTOFGeometry::NStripB();
1159   Int_t nStripC = AliTOFGeometry::NStripC();
1160
1161   // Int_t isector = in[0];
1162   Int_t iplate = in[1];
1163   Int_t istrip = in[2];
1164   //Int_t ipadX = in[3]; 
1165   //Int_t ipadZ = in[4]; 
1166   
1167   Int_t stripOffset = 0;
1168   switch (iplate) {
1169   case 0:
1170     stripOffset = 0;
1171       break;
1172   case 1:
1173     stripOffset = nStripC;
1174     break;
1175   case 2:
1176     stripOffset = nStripC+nStripB;
1177     break;
1178   case 3:
1179     stripOffset = nStripC+nStripB+nStripA;
1180     break;
1181   case 4:
1182     stripOffset = nStripC+nStripB+nStripA+nStripB;
1183     break;
1184   default:
1185     //AliDebug(2,Form("Wrong plate number in TOF (%d) !",iplate));
1186       stripOffset=-1;
1187       break;
1188   };
1189   
1190   if (stripOffset<0 || stripOffset>92) return -1;
1191   else 
1192       return (stripOffset+istrip);
1193   
1194 }
1195 //---------------------------------------------------------------
1196 Bool_t  AliTOFQADataMakerRec::CheckVolumeID(const Int_t * const volumeID)
1197 {
1198     //
1199     //Checks volume ID validity
1200     //   
1201     for (Int_t j=0;j<5;j++){
1202         if (volumeID[j]<0) {
1203           //AliDebug(2,Form("Invalid detector volume index for volumeID[%i]",j));
1204             return kFALSE;
1205         }
1206     }
1207     return kTRUE; 
1208 }
1209
1210 //---------------------------------------------------------------
1211 Bool_t  AliTOFQADataMakerRec::CheckEquipID(const Int_t * const equipmentID)
1212 {
1213     //
1214     //Checks equipment ID validity    
1215    for (Int_t j=0;j<5;j++){
1216         if (equipmentID[j]<0) {
1217           // AliDebug(2,Form("Invalid equipment volume index for equipmentID[%i]",j));
1218           return kFALSE;
1219         }
1220    }
1221    return kTRUE;
1222 }
1223 //---------------------------------------------------------------
1224 Bool_t  AliTOFQADataMakerRec::FilterLTMData(const Int_t * const equipmentID) const
1225 {
1226   /*It returns kTRUE if data come from LTM.
1227     It thus filters trigger-related signals  */
1228
1229   Int_t ddl, trm, tdc;
1230   //if (!CheckEquipID(equipmentID)) return kFALSE;
1231   ddl = equipmentID[0];
1232   trm = equipmentID[1];
1233   tdc = equipmentID[3];
1234   
1235   if ((ddl%2==1) && (trm==3) && (tdc>11 && tdc<15))
1236     return kTRUE;
1237   else 
1238     return kFALSE;
1239  
1240 }
1241 //---------------------------------------------------------------
1242 Bool_t  AliTOFQADataMakerRec::FilterSpare(const Int_t * const equipmentID) const
1243 {
1244   /*It returns kTRUE if data come from spare 
1245     equipment ID. 
1246     So far only check on TRM 3 crate left is implemented */
1247
1248   Int_t ddl, trm, tdc;
1249   //if (!CheckEquipID(equipmentID)) return kFALSE;
1250   ddl = equipmentID[0];
1251   trm = equipmentID[1];
1252   tdc = equipmentID[3];
1253   
1254   if ((ddl%2==1) && (trm==3) && (tdc>2 && tdc<12))
1255     return kTRUE;
1256   else 
1257     return kFALSE;
1258 }
1259
1260 //-----------------------------------------------------------------------------
1261 void AliTOFQADataMakerRec::GetGeo2LTMIndex(const Int_t * const detind, Int_t *indexLTM) {
1262   //
1263   // getting LTMmatrix indexes for current digit
1264   //
1265   Int_t stripId=GetStripIndex(detind);
1266
1267   if (detind[1]==0 || detind[1]==1 || (detind[1]==2 && detind[2]<=7)) { //A side
1268     if (detind[4]<24){ //R
1269       indexLTM[0] = detind[0]*2;
1270     } else { //L
1271       indexLTM[0] = detind[0]*2+1;
1272     }  
1273     indexLTM[1]=stripId;
1274
1275   } else { //C side
1276     if (detind[4]<24){
1277       indexLTM[0] = detind[0]*2+36;
1278     } else {
1279       indexLTM[0] = (detind[0]*2+1)+36;
1280     }
1281     indexLTM[1]=90-stripId; 
1282   }
1283   
1284   // if (indexLTM[0]<36) { //A side
1285   //   if (detind[1] ==0){
1286   //     indexLTM[1] = detind[2];
1287   //   }
1288   //   else if (detind[1] ==1){
1289   //     indexLTM[1] = detind[2]+nStripB;
1290   //   }
1291   //   else if (detind[1] ==2){
1292   //     indexLTM[1] = detind[2]+19*2;
1293   //   }
1294   //   else{
1295   //     AliError("Smth Wrong!!!");
1296   //   }
1297   // }
1298   // else { //C side
1299   //   if (detind[1]==2){
1300   //     if (detind[4]<24)
1301   //    indexLTM[1] = (nStripAL-detind[2])+nStripC+nStripB;
1302   //     else 
1303   //    indexLTM[1] = (nStripAR-detind[2])+nStripC+nStripB;
1304   //   }
1305   //   else if (detind[1] ==3){
1306   //     indexLTM[1] = (nStripB-detind[2])+nStripC;
1307   //   }
1308   //   else if (detind[1] ==4){
1309   //     indexLTM[1] = nStripC-detind[2];
1310   //   }
1311   //   else{
1312   //     AliError("Smth Wrong!!!");
1313   //   }
1314   // }  
1315 }
1316
1317 //-----------------------------------------------------------------------------
1318 void AliTOFQADataMakerRec::GetGeo2CTTMIndex(Int_t *detind, Int_t *indexCTTM) {
1319   //
1320   // Returns CTTM index corresponding to the detector element detind
1321   //
1322   GetGeo2LTMIndex(detind,indexCTTM);
1323   indexCTTM[1]/=2;
1324   return;
1325 }
1326
1327 //-------------------------------------------------------------------------
1328 Int_t AliTOFQADataMakerRec::GetNumberOfFiredMacropad(AliRawReader * rawReader){
1329   
1330   Int_t nFired=0;
1331   TClonesArray * clonesRawData;  
1332   if (!rawReader) return 0;  
1333   for (Int_t iDDL = 0; iDDL < AliTOFGeometry::NDDL()*AliTOFGeometry::NSectors(); iDDL++){
1334     rawReader->Reset();
1335     fTOFRawStream.LoadRawDataBuffersV2(iDDL); 
1336     clonesRawData = (TClonesArray*)fTOFRawStream.GetRawData();
1337     for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
1338       AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);    
1339       if (tofRawDatum->GetTOF()){       
1340         if ( (tofRawDatum->GetTRM()==3)&&
1341              (tofRawDatum->GetTDC()>11)&&
1342              (tofRawDatum->GetTDC()<15)) {        
1343           nFired+=1;  
1344         }
1345       }
1346     }    
1347   }//loop over DDLs
1348   return nFired; 
1349 }
1350
1351 //----------------------------------------------------------------
1352 void AliTOFQADataMakerRec::GetCTTMIndex(Int_t *equipid, Int_t *indexCTTM) {
1353   //
1354   // Returns CTTM index corresponding to the equipment id equipid, only for LTM hits
1355   // equipid = (crate, trm, chain, tdc, channel)
1356
1357   if ( (equipid[1]!=3)||(equipid[3]<12) ){
1358     indexCTTM[0]=-1;
1359     indexCTTM[1]=-1;
1360     return;
1361   }  
1362   Int_t modStrip2LTM[3][8]={ { 0, 1, 2, 3, 4, 5, 6, 7},
1363                              { 8, 9, 10, 11, 12, 13, 14, 15},
1364                              {16, 17, 18, 19, 20, 21, 22, 23}
1365                             }; 
1366
1367   Int_t DDL2LTMmatrix[72]={0,1,37,36,2,3,39,38,4,5,41,40,6,7,43,42,8,9,45,44,10,11,47,46,12,13,49,48,14,15,51,50,16,17,53,52,18,19,
1368                            55,54,20,21,57,56,22,23,59,58,24,25,61,60,26,27,63,62,28,29,65,64,30,31,67,66,32,33,69,68,34,35,71,70};
1369
1370   Int_t itdc=equipid[3]%12;
1371   Int_t crate=-1;
1372   if (equipid[2]==0)
1373    crate=equipid[0]-1;
1374   else crate=equipid[0];
1375   
1376   indexCTTM[0]=DDL2LTMmatrix[crate];
1377   indexCTTM[1]=modStrip2LTM[itdc][equipid[4]];      
1378   return;
1379 }
1380
1381 //_____________________________________________________________________________
1382 void AliTOFQADataMakerRec::ReadHistogramRangeFromFile(const Char_t * filename)
1383 {
1384   //
1385   // read histogram ranges from configuration file
1386   //
1387   if (!filename) {
1388     AliInfo("Config file with histograms ranges not found or invalid -> use default values.");
1389     SetDefaultHistogramRange();
1390     SetDefaultCutNmaxFiredMacropad();
1391     return;
1392   }
1393   
1394   std::fstream configFile;
1395   configFile.open(filename, std::fstream::in);
1396   if (!configFile.is_open()){
1397     AliInfo("Cannot open config file with histograms ranges -> use default values.");
1398     SetDefaultHistogramRange();
1399     return;
1400   }
1401   
1402   //check file size
1403   Int_t begin = configFile.tellg();
1404   configFile.seekg(0, std::fstream::end); /* end */
1405   Int_t end = configFile.tellg();
1406   Int_t size = end - begin;
1407   configFile.seekg(0, std::fstream::beg); /* rewind file */
1408   if (size <= 0){
1409     AliInfo(Form("Unexpected EOF of config file with histograms ranges. File size: %d -> use default values", size));
1410     SetDefaultHistogramRange();
1411     return;
1412   }
1413   
1414   Int_t minMulti=9999, maxMulti=-9999;
1415   Int_t nbinsMulti=0,nbinsTime=0;
1416   Float_t minTime=9999.0, maxTime=-9999.0;
1417   Int_t cutFiredMacropad=0;
1418   TString endoflist;
1419   while (!configFile.eof()) {
1420     configFile >> cutFiredMacropad >> minMulti >> maxMulti >> minTime >> maxTime;
1421     configFile >> endoflist;
1422     if (endoflist.Contains("end")) break;
1423   }
1424
1425   //set multiplicity histo ranges
1426   if (minMulti>maxMulti){
1427     AliInfo("Invalid range for multiplicity histogram set. Changing to default values.");
1428     SetDefaultMultiHistogramRange();
1429   } else {
1430     nbinsMulti = maxMulti-minMulti;
1431     SetNbinsMultiplicityHisto(nbinsMulti);
1432     SetMultiplicityHistoRange(minMulti,maxMulti);
1433     //AliInfo(Form("Setting multiplicity histogram ranges to: multMin = %i - multMax = %i - nMultBins = %i", fgRangeMinMultiplicity, fgRangeMaxMultiplicity, fgNbinsMultiplicity));
1434   }
1435
1436   //set time histo ranges
1437   if (minTime>maxTime){
1438     AliInfo("Invalid range for time histogram set. Changing to default values.");
1439     SetDefaultTimeHistogramRange();
1440   } else {
1441     nbinsTime = TMath::Nint((maxTime - minTime)/fgkNbinsWidthTime);//ns
1442     maxTime=minTime+nbinsTime*fgkNbinsWidthTime;//ns
1443     SetNbinsTimeHisto(nbinsTime);
1444     SetTimeHistoRange(minTime,maxTime);
1445     //AliInfo(Form("Setting time histogram ranges to: timeMin = %5.2f ns - timeMax = %5.2f ns - nTimeBins = %i", fgRangeMinTime, fgRangeMaxTime,fgNbinsTime));
1446   } 
1447  
1448   if ((cutFiredMacropad>0)&&(cutFiredMacropad<fgkFiredMacropadLimit)){
1449     AliInfo("Invalid value for cut on fired macropad. Changing to default values.");
1450     SetDefaultCutNmaxFiredMacropad();
1451   } else {
1452     SetCutNmaxFiredMacropad(cutFiredMacropad);
1453     //AliInfo(Form("Setting cut on fired macropad to:  = %i",cutFiredMacropad));
1454   } 
1455   AliInfo(Form("Setting: multMin = %i - multMax = %i - nMultBins = %i, timeMin = %5.2f ns - timeMax = %5.2f ns - nTimeBins = %i, cutMaxFiredMacropad = %i", 
1456                fgRangeMinMultiplicity, fgRangeMaxMultiplicity, fgNbinsMultiplicity, fgRangeMinTime, fgRangeMaxTime,fgNbinsTime, cutFiredMacropad));
1457   configFile.close();
1458   return;
1459 }
1460
1461 //_____________________________________________________________________________
1462 void AliTOFQADataMakerRec::SetDefaultHistogramRange()
1463 {
1464   //
1465   // set default histogram ranges (tuned on 2011 pp collisions)
1466   // 
1467   //AliInfo("Setting all histogram ranges to default values.");
1468   SetDefaultMultiHistogramRange();
1469   SetDefaultTimeHistogramRange();
1470   SetDefaultCutNmaxFiredMacropad();
1471   return;
1472 }
1473
1474 //_____________________________________________________________________________
1475 void AliTOFQADataMakerRec::SetDefaultMultiHistogramRange()
1476 {
1477   //
1478   // set default histogram ranges (tuned on 2011 pp collisions)
1479   // 
1480   SetMultiplicityHistoRange (0, 200);
1481   SetNbinsMultiplicityHisto(200);
1482   //AliInfo("Setting Multiplicity histogram ranges to default values.");
1483   //AliInfo(Form("multMin = %i - multMax = %i - nMultBins = %i",
1484   //           fgRangeMinMultiplicity, fgRangeMaxMultiplicity, fgNbinsMultiplicity));
1485   return;
1486 }
1487
1488 //_____________________________________________________________________________
1489 void AliTOFQADataMakerRec::SetDefaultTimeHistogramRange()
1490 {
1491   //
1492   // set default histogram ranges (tuned on 2011 pp collisions)
1493   // 
1494   SetNbinsTimeHisto(250);
1495   SetTimeHistoRange (0.0,610.);   
1496   
1497   // AliInfo("Setting Time histogram ranges to default values:");
1498   // AliInfo(Form("timeMin = %5.2f ns - timeMax = %5.2f ns - nTimeBins = %i",
1499   //           fgRangeMinTime, fgRangeMaxTime,fgNbinsTime));
1500   return;
1501 }
1502
1503 //------------------------------------------------------------------
1504 void AliTOFQADataMakerRec::SetDefaultCutNmaxFiredMacropad()
1505 {
1506   //
1507   // set default cut on fired macropad 
1508   // 
1509   SetCutNmaxFiredMacropad(50); 
1510   // AliInfo(Form("Setting cut on fired macropad to default values: NfiredMacropad = %i", fgCutNmaxFiredMacropad));
1511   return;
1512 }