]> git.uio.no Git - u/mrichter/AliRoot.git/blob - VZERO/AliVZEROQADataMakerRec.cxx
Correction of additional compilation errors from the first compilation trial.
[u/mrichter/AliRoot.git] / VZERO / AliVZEROQADataMakerRec.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 //  Produces the data needed to calculate the quality assurance 
18 //  All data must be mergeable objects
19 //  Handles ESDs and Raws
20 //  Histos defined will be used for Raw Data control and monitoring
21
22 // --- ROOT system ---
23 #include <TClonesArray.h>
24 #include <TFile.h> 
25 #include <TF1.h> 
26 #include <TH1F.h> 
27 #include <TH1I.h> 
28 #include <TH2I.h> 
29 #include <TH2F.h> 
30 #include <TGraph.h> 
31 #include <TParameter.h>
32 #include <TTimeStamp.h>
33
34 // --- Standard library ---
35
36 // --- AliRoot header files ---
37 #include "AliESDEvent.h"
38 #include "AliLog.h"
39 #include "AliCDBManager.h"
40 #include "AliCDBStorage.h"
41 #include "AliCDBEntry.h"
42 #include "AliVZEROQADataMakerRec.h"
43 #include "AliQAChecker.h"
44 #include "AliRawReader.h"
45 #include "AliVZERORawStream.h"
46 #include "AliVZEROdigit.h"
47 #include "AliVZEROConst.h"
48 #include "AliVZEROReconstructor.h"
49 #include "AliVZEROTrending.h"
50 #include "AliVZEROCalibData.h"
51 #include "AliCTPTimeParams.h"
52 #include "event.h"
53
54 const Float_t kMinBBA = 68. ;
55 const Float_t kMaxBBA = 100. ;
56 const Float_t kMinBBC = 75.5 ;
57 const Float_t kMaxBBC = 100. ;
58 const Float_t kMinBGA = 54. ;
59 const Float_t kMaxBGA = 58. ;
60 const Float_t kMinBGC = 69.5 ;
61 const Float_t kMaxBGC = 74. ;
62
63  
64  
65  
66
67 ClassImp(AliVZEROQADataMakerRec)
68            
69 //____________________________________________________________________________ 
70 AliVZEROQADataMakerRec::AliVZEROQADataMakerRec() : 
71 AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kVZERO), "VZERO Quality Assurance Data Maker"),
72   fCalibData(0x0),
73 //  fEvent(0), 
74 //  fNTotEvents(0), 
75 //  fNSubEvents(0), 
76 //  fTrendingUpdateEvent(0), 
77 //  fNTrendingUpdates(0), 
78   fTrendingUpdateTime(0), 
79   fCycleStartTime(0), 
80   fCycleStopTime(0),
81   fTimeSlewing(0)
82     
83 {
84   // Constructor
85    
86   AliDebug(AliQAv1::GetQADebugLevel(), "Construct VZERO QA Object");
87
88   for(Int_t i=0; i<64; i++){  
89     fEven[i] = 0;   
90     fOdd[i]  = 0;
91   }
92   
93   for(Int_t i=0; i<128; i++){  
94     fADCmean[i] = 0.0;   }      
95 }
96
97 //____________________________________________________________________________ 
98 AliVZEROQADataMakerRec::AliVZEROQADataMakerRec(const AliVZEROQADataMakerRec& qadm) :
99   AliQADataMakerRec(),
100   fCalibData(0x0),
101   //  fEvent(0), 
102   //  fNTotEvents(0), 
103   //  fNSubEvents(0), 
104   //  fTrendingUpdateEvent(0), 
105   //  fNTrendingUpdates(0), 
106   fTrendingUpdateTime(0), 
107   fCycleStartTime(0), 
108   fCycleStopTime(0),
109   fTimeSlewing(0)
110   
111 {
112   // Copy constructor 
113   
114   SetName((const char*)qadm.GetName()) ; 
115   SetTitle((const char*)qadm.GetTitle()); 
116 }
117
118 //__________________________________________________________________
119 AliVZEROQADataMakerRec& AliVZEROQADataMakerRec::operator = (const AliVZEROQADataMakerRec& qadm )
120 {
121   // Equal operator
122   
123   this->~AliVZEROQADataMakerRec();
124   new(this) AliVZEROQADataMakerRec(qadm);
125   return *this;
126 }
127
128 //____________________________________________________________________________
129 AliVZEROCalibData* AliVZEROQADataMakerRec::GetCalibData() const
130
131 {
132   AliCDBManager *man = AliCDBManager::Instance();
133
134   AliCDBEntry *entry=0;
135
136   entry = man->Get("VZERO/Calib/Data",fRun);
137   if(!entry){
138     AliWarning("Load of calibration data from default storage failed!");
139     AliWarning("Calibration data will be loaded from local storage ($ALICE_ROOT)");
140         
141     man->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
142     entry = man->Get("VZERO/Calib/Data",fRun);
143   }
144   // Retrieval of data in directory VZERO/Calib/Data:
145
146   AliVZEROCalibData *calibdata = 0;
147
148   if (entry) calibdata = (AliVZEROCalibData*) entry->GetObject();
149   if (!calibdata)  AliFatal("No calibration data from calibration database !");
150
151   return calibdata;
152 }
153
154 //____________________________________________________________________________ 
155 void AliVZEROQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
156 {
157   // Detector specific actions at end of cycle
158   // Does the QA checking
159   ResetEventTrigClasses();
160   //
161   AliQAChecker::Instance()->Run(AliQAv1::kVZERO, task, list) ;
162   
163   if(task == AliQAv1::kRAWS){
164     TTimeStamp currentTime;
165     fCycleStopTime = currentTime.GetSec();
166   }
167
168   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
169     if (! IsValidEventSpecie(specie, list)) continue ;
170     SetEventSpecie(AliRecoParam::ConvertIndex(specie));
171     if(task == AliQAv1::kRAWS) {
172     } else if (task == AliQAv1::kESDS) {
173     }
174   }
175 }
176
177 //____________________________________________________________________________ 
178 void AliVZEROQADataMakerRec::InitESDs()
179 {
180   // Creates histograms to control ESDs
181   
182   const Bool_t expert   = kTRUE ; 
183   const Bool_t image    = kTRUE ; 
184         
185   TH2F * h2d;
186   TH1I * h1i;
187   TH1F * h1d;
188                 
189   h1i = new TH1I("H1I_Cell_Multiplicity_V0A", "Cell Multiplicity in V0A;Multiplicity (Nb of Cell);Counts", 35, 0, 35) ;  
190   Add2ESDsList(h1i, kCellMultiV0A, !expert, image)  ;  
191                                                                                                         
192   h1i = new TH1I("H1I_Cell_Multiplicity_V0C", "Cell Multiplicity in V0;Multiplicity (Nb of Cell);Counts", 35, 0, 35) ;  
193   Add2ESDsList(h1i, kCellMultiV0C, !expert, image)  ;  
194    
195   h1d = new TH1F("H1D_MIP_Multiplicity_V0A", "MIP Multiplicity in V0A;Multiplicity (Nb of MIP);Counts", 1000, 0, 1000) ;  
196   Add2ESDsList(h1d, kMIPMultiV0A, !expert, image)  ;  
197   
198   h1d = new TH1F("H1D_MIP_Multiplicity_V0C", "MIP Multiplicity in V0C;Multiplicity (Nb of MIP);Counts", 1000, 0, 1000) ;  
199   Add2ESDsList(h1d, kMIPMultiV0C, !expert, image)  ;  
200
201   h2d = new TH2F("H2D_MIP_Multiplicity_Channel", "MIP Multiplicity per Channel;Channel;Multiplicity (Nb of MIP)",64, 0, 64, 100, 0, 100) ;  
202   Add2ESDsList(h2d, kMIPMultiChannel, !expert, image)  ;  
203   
204   h1d = new TH1F("H1D_BBFlag_Counters", "BB Flag Counters;Channel;Counts",64, 0, 64) ;  
205   Add2ESDsList(h1d, kBBFlag, !expert, image)  ;  
206   
207   h1d = new TH1F("H1D_BGFlag_Counters", "BG Flag Counters;Channel;Counts",64, 0, 64) ;  
208   Add2ESDsList(h1d, kBGFlag, !expert, image)  ;  
209   
210   h2d = new TH2F("H2D_Charge_Channel", "ADC Charge per channel;Channel;Charge (ADC counts)",64, 0, 64, 1024, 0, 1024) ;  
211   Add2ESDsList(h2d, kChargeChannel, !expert, image)  ;  
212   
213   h2d = new TH2F("H2D_Time_Channel", "Time per channel;Channel;Time (ns)",64, 0, 64, 400, -100, 100) ;  
214   Add2ESDsList(h2d, kTimeChannel, !expert, image)  ;  
215   
216   h1d = new TH1F("H1D_V0A_Time", "Mean V0A Time;Time (ns);Counts",1000, -100., 100.);
217   Add2ESDsList(h1d,kESDV0ATime, !expert, image); 
218   
219   h1d = new TH1F("H1D_V0C_Time", "Mean V0C Time;Time (ns);Counts",1000, -100., 100.);
220   Add2ESDsList(h1d,kESDV0CTime, !expert, image); 
221   
222   h1d = new TH1F("H1D_Diff_Time", "Diff Time V0A - V0C;Diff Time V0A - V0C (ns);Counts",1000, -200., 200.);
223   Add2ESDsList(h1d,kESDDiffTime, !expert, image); 
224   //
225   ClonePerTrigClass(AliQAv1::kESDS); // this should be the last line    
226 }
227
228 //____________________________________________________________________________ 
229 void AliVZEROQADataMakerRec::InitRaws()
230 {
231   // Creates RAW histograms in Raws subdir
232
233   const Bool_t expert   = kTRUE ; 
234   const Bool_t saveCorr = kTRUE ; 
235   const Bool_t image    = kTRUE ; 
236
237   const Int_t kNintegrator  =    2;
238  
239   const Int_t kNTdcTimeBins  = 1280;
240   const Float_t kTdcTimeMin    =    0.;
241   const Float_t kTdcTimeMax    = 75.;
242     const Int_t kNTdcWidthBins =  128;
243   const Float_t kTdcWidthMin   =    0;
244   const Float_t kTdcWidthMax   =  50.;
245   const Int_t kNChargeBins   = 1024;
246   const Float_t kChargeMin     =    0;
247   const Float_t kChargeMax     = 1024;
248   const Int_t kNChannelBins  =   64;
249   const Float_t kChannelMin    =    0;
250   const Float_t kChannelMax    =   64;
251   const Int_t kNPedestalBins =  200;
252   const Float_t kPedestalMin   =    0;
253   const Float_t kPedestalMax   =  200;
254   const Int_t kNMIPBins      = 512;
255   const Float_t kMIPMin        =   0;
256   const Float_t kMIPMax        = 16;
257
258   TH2I * h2i;
259   TH2F * h2d;
260   TH1I * h1i;
261   TH1F * h1d;
262
263   int iHisto =0;
264   // Creation of Trigger Histogram
265   h1d = new TH1F("H1D_Trigger_Type", "V0 Trigger Type;;Counts", 4,0 ,4) ;  
266   Add2RawsList(h1d,kTriggers, !expert, image, saveCorr);   iHisto++;
267         h1d->SetFillColor(29);
268         h1d->SetLineWidth(2);
269         h1d->GetXaxis()->SetLabelSize(0.06);
270         h1d->GetXaxis()->SetNdivisions(808,kFALSE);
271         h1d->GetXaxis()->SetBinLabel(1, "V0-AND");
272         h1d->GetXaxis()->SetBinLabel(2, "V0-OR");
273         h1d->GetXaxis()->SetBinLabel(3, "V0-BGA");
274         h1d->GetXaxis()->SetBinLabel(4, "V0-BGC");
275
276         h2d = new TH2F("H2D_Trigger_Type", "V0 Trigger Type;V0A;V0C", 4,0 ,4,4,0,4) ;  
277         Add2RawsList(h2d,kTriggers2, !expert, image, saveCorr);   iHisto++;
278         h2d->SetDrawOption("coltext");
279         h2d->GetXaxis()->SetLabelSize(0.06);
280         h2d->GetXaxis()->SetNdivisions(808,kFALSE);
281         h2d->GetYaxis()->SetNdivisions(808,kFALSE);
282         h2d->GetXaxis()->SetBinLabel(1, "Empty");
283         h2d->GetXaxis()->SetBinLabel(2, "Fake");
284         h2d->GetXaxis()->SetBinLabel(3, "BB");
285         h2d->GetXaxis()->SetBinLabel(4, "BG");
286         h2d->GetYaxis()->SetBinLabel(1, "Empty");
287         h2d->GetYaxis()->SetBinLabel(2, "Fake");
288         h2d->GetYaxis()->SetBinLabel(3, "BB");
289         h2d->GetYaxis()->SetBinLabel(4, "BG");
290
291   // Creation of Cell Multiplicity Histograms
292   h1i = new TH1I("H1I_Multiplicity_V0A", "Cell Multiplicity in V0A;# of Cells;Entries", 35, 0, 35) ;  
293   Add2RawsList(h1i,kMultiV0A, expert, image, saveCorr);   iHisto++;
294   h1i = new TH1I("H1I_Multiplicity_V0C", "Cell Multiplicity in V0C;# of Cells;Entries", 35, 0, 35) ;  
295   Add2RawsList(h1i,kMultiV0C, expert, image, saveCorr);   iHisto++;
296  
297   // Creation of Total Charge Histograms
298   h1d = new TH1F("H1D_Charge_V0A", "Total Charge in V0A;Charge [ADC counts];Counts", 2000, 0, 10000) ;  
299   Add2RawsList(h1d,kChargeV0A, expert, !image, saveCorr);   iHisto++;
300   h1d = new TH1F("H1D_Charge_V0C", "Total Charge in V0C;Charge [ADC counts];Counts", 2000, 0, 10000) ;  
301   Add2RawsList(h1d,kChargeV0C, expert, !image, saveCorr);   iHisto++;
302   h1d = new TH1F("H1D_Charge_V0", "Total Charge in V0;Charge [ADC counts];Counts", 2000, 0, 20000) ;  
303   Add2RawsList(h1d,kChargeV0, expert, !image, saveCorr);   iHisto++;
304   
305   // Creation of MIP Histograms
306   h1d = new TH1F("H1D_MIP_V0A", "Total MIP in V0A;Multiplicity [MIP];Counts", kNMIPBins,kMIPMin ,32*kMIPMax) ;  
307   Add2RawsList(h1d,kRawMIPV0A, expert, !image, saveCorr);   iHisto++;
308   h1d = new TH1F("H1D_MIP_V0C", "Total MIP in V0C;Multiplicity [MIP];Counts", kNMIPBins,kMIPMin ,32*kMIPMax) ;  
309   Add2RawsList(h1d,kRawMIPV0C, expert, !image, saveCorr);   iHisto++;
310   h1d = new TH1F("H1D_MIP_V0", "Total MIP in V0;Multiplicity [MIP];Counts", 2*kNMIPBins,kMIPMin ,64*kMIPMax) ;  
311   Add2RawsList(h1d,kRawMIPV0, expert, !image, saveCorr);   iHisto++;
312   h2d = new TH2F("H2D_MIP_Channel", "Nb of MIP per channel;Channel;# of Mips", kNChannelBins, kChannelMin, kChannelMax,kNMIPBins,kMIPMin ,kMIPMax) ;  
313   Add2RawsList(h2d,kRawMIPChannel, expert, !image, !saveCorr);   iHisto++;
314   
315  
316
317   // Creation of Charge EoI histogram 
318   h2d = new TH2F("H2D_ChargeEoI", "Charge Event of Interest;Channel Number;Charge [ADC counts]"
319                  ,kNChannelBins, kChannelMin, kChannelMax, kNChargeBins, kChargeMin, kChargeMax);
320   Add2RawsList(h2d,kChargeEoI, !expert, image, !saveCorr); iHisto++;
321
322   for(Int_t iInt=0;iInt<kNintegrator;iInt++){
323     // Creation of Pedestal histograms 
324     h2i = new TH2I(Form("H2I_Pedestal_Int%d",iInt), Form("Pedestal (Int%d);Channel;Pedestal [ADC counts]",iInt)
325                    ,kNChannelBins, kChannelMin, kChannelMax,kNPedestalBins,kPedestalMin ,kPedestalMax );
326     Add2RawsList(h2i,(iInt == 0 ? kPedestalInt0 : kPedestalInt1), expert, !image, !saveCorr); iHisto++;
327         
328
329     // Creation of Charge EoI histograms 
330     h2i = new TH2I(Form("H2I_ChargeEoI_Int%d",iInt), Form("Charge EoI (Int%d);Channel;Charge [ADC counts]",iInt)
331                    ,kNChannelBins, kChannelMin, kChannelMax, kNChargeBins, kChargeMin, kChargeMax);
332     Add2RawsList(h2i,(iInt == 0 ? kChargeEoIInt0 : kChargeEoIInt1), expert, image, !saveCorr); iHisto++;
333     
334     h2i = new TH2I(Form("H2I_ChargeEoI_BB_Int%d",iInt), Form("Charge EoI w/ BB Flag (Int%d);Channel;Charge [ADC counts]",iInt)
335                    ,kNChannelBins, kChannelMin, kChannelMax, kNChargeBins, kChargeMin, kChargeMax);
336     Add2RawsList(h2i,(iInt == 0 ? kChargeEoIBBInt0 : kChargeEoIBBInt1), expert, !image, !saveCorr); iHisto++;
337     
338     h2i = new TH2I(Form("H2I_ChargeEoI_BG_Int%d",iInt), Form("Charge EoI w/ BG Flag (Int%d);Channel;Charge [ADC counts]",iInt)
339                    ,kNChannelBins, kChannelMin, kChannelMax, kNChargeBins, kChargeMin, kChargeMax);
340     Add2RawsList(h2i,(iInt == 0 ?  kChargeEoIBGInt0: kChargeEoIBGInt1), expert, !image, !saveCorr); iHisto++;
341
342     // Creation of Charge versus LHC Clock histograms 
343     h2d = new TH2F(Form("H2D_ChargeVsClock_Int%d",iInt), Form("Charge Versus LHC-Clock (Int%d);Channel;LHCClock;Charge [ADC counts]",iInt)
344                    ,kNChannelBins, kChannelMin, kChannelMax,21, -10.5, 10.5 );
345     Add2RawsList(h2d,(iInt == 0 ? kChargeVsClockInt0 : kChargeVsClockInt1 ), expert, !image, !saveCorr); iHisto++;
346         
347     // Creation of Minimum Bias Charge histograms 
348     for(Int_t iBB=0;iBB<2;iBB++){
349       for(Int_t iBG=0;iBG<2;iBG++){
350         h2i = new TH2I(Form("H2I_ChargeMB_BB%d_BG%d_Int%d",iBB,iBG,iInt), Form("MB Charge (BB=%d, BG=%d, Int=%d);Channel;Charge [ADC counts]",iBB,iBG,iInt)
351                        ,kNChannelBins, kChannelMin, kChannelMax,kNChargeBins, kChargeMin, kChargeMax);
352         int idx;
353         if(iInt==0){
354           if(iBB==0){
355             if(iBG==0) idx = kChargeMBBB0BG0Int0;
356             else idx = kChargeMBBB0BG1Int0;
357           } else {
358             if(iBG==0) idx = kChargeMBBB1BG0Int0;
359             else idx = kChargeMBBB1BG1Int0;
360           }
361         } else {
362           if(iBB==0){
363             if(iBG==0) idx = kChargeMBBB0BG0Int1;
364             else idx = kChargeMBBB0BG1Int1;
365           } else {
366             if(iBG==0) idx = kChargeMBBB1BG0Int1;
367             else idx = kChargeMBBB1BG1Int1;
368           }
369         }
370         Add2RawsList(h2i,idx, expert, !image, !saveCorr); iHisto++;
371       }
372     }
373         
374   }
375  
376   // Creation of Time histograms 
377   h2i = new TH2I("H2I_Width", "HPTDC Width;Channel;Width [ns]",kNChannelBins, kChannelMin, kChannelMax, kNTdcWidthBins, kTdcWidthMin, kTdcWidthMax);
378   Add2RawsList(h2i,kWidth, expert, !image, !saveCorr); iHisto++;
379
380   h2i = new TH2I("H2I_Width_BB", "HPTDC Width w/ BB Flag condition;Channel;Width [ns]",kNChannelBins, kChannelMin, kChannelMax, kNTdcWidthBins, kTdcWidthMin, kTdcWidthMax);
381   Add2RawsList(h2i,kWidthBB, expert, !image, !saveCorr); iHisto++;
382
383   h2i = new TH2I("H2I_Width_BG", "HPTDC Width w/ BG Flag condition;Channel;Width [ns]",kNChannelBins, kChannelMin, kChannelMax, kNTdcWidthBins, kTdcWidthMin, kTdcWidthMax);
384   Add2RawsList(h2i,kWidthBG, expert, !image, !saveCorr); iHisto++;
385
386   h2i = new TH2I("H2I_HPTDCTime", "HPTDC Time;Channel;Leading Time [ns]",kNChannelBins, kChannelMin, kChannelMax, kNTdcTimeBins, kTdcTimeMin, kTdcTimeMax);
387   Add2RawsList(h2i,kHPTDCTime, expert, image, !saveCorr); iHisto++;
388
389   h2i = new TH2I("H2I_HPTDCTime_BB", "HPTDC Time w/ BB Flag condition;Channel;Leading Time [ns]",kNChannelBins, kChannelMin, kChannelMax, kNTdcTimeBins, kTdcTimeMin, kTdcTimeMax);
390   Add2RawsList(h2i,kHPTDCTimeBB, !expert, image, !saveCorr); iHisto++;
391
392   h2i = new TH2I("H2I_HPTDCTime_BG", "HPTDC Time w/ BG Flag condition;Channel;Leading Time [ns]",kNChannelBins, kChannelMin, kChannelMax, kNTdcTimeBins, kTdcTimeMin, kTdcTimeMax);
393   Add2RawsList(h2i,kHPTDCTimeBG, !expert, image, !saveCorr); iHisto++;
394         
395   h1d = new TH1F("H1D_V0A_Time", "V0A Time;Time [ns];Counts",kNTdcTimeBins, kTdcTimeMin, kTdcTimeMax);
396   Add2RawsList(h1d,kV0ATime, expert, !image, saveCorr); iHisto++;
397         
398   h1d = new TH1F("H1D_V0C_Time", "V0C Time;Time [ns];Counts",kNTdcTimeBins, kTdcTimeMin, kTdcTimeMax);
399   Add2RawsList(h1d,kV0CTime, expert, !image, saveCorr); iHisto++;
400         
401   h1d = new TH1F("H1D_Diff_Time","Diff V0A-V0C Time;Time [ns];Counts",kNTdcTimeBins, -50., 50.);
402   Add2RawsList(h1d,kDiffTime, expert, image, saveCorr); iHisto++;
403
404     h2d = new TH2F("H2D_TimeV0A_V0C", "Mean Time in V0C versus V0A;Time V0A [ns];Time V0C [ns]", 
405                 kNTdcTimeBins/8, kTdcTimeMin,kTdcTimeMax,kNTdcTimeBins/8, kTdcTimeMin,kTdcTimeMax) ;  
406     Add2RawsList(h2d,kTimeV0AV0C, !expert, image, !saveCorr);   iHisto++;
407         
408   // Creation of Flag versus LHC Clock histograms 
409
410   h1d = new TH1F("H1D_BBFlagPerChannel", "BB-Flags Versus Channel;Channel;BB Flags Count",kNChannelBins, kChannelMin, kChannelMax );
411   h1d->SetMinimum(0);
412   Add2RawsList(h1d,kBBFlagsPerChannel, !expert, image, !saveCorr); iHisto++;
413
414   h2d = new TH2F("H2D_BBFlagVsClock", "BB-Flags Versus LHC-Clock;Channel;LHC Clocks",kNChannelBins, kChannelMin, kChannelMax,21, -10.5, 10.5 );
415   Add2RawsList(h2d,kBBFlagVsClock, expert, !image, !saveCorr); iHisto++;
416         
417   h2d = new TH2F("H2D_BGFlagVsClock", "BG-Flags Versus LHC-Clock;Channel;LHC Clocks",kNChannelBins, kChannelMin, kChannelMax,21, -10.5, 10.5 );
418   Add2RawsList(h2d,kBGFlagVsClock, expert, !image, !saveCorr); iHisto++;
419   // 
420   AliDebug(AliQAv1::GetQADebugLevel(), Form("%d Histograms has been added to the Raws List",iHisto));
421   //
422   ClonePerTrigClass(AliQAv1::kRAWS); // this should be the last line
423 }
424
425 //____________________________________________________________________________ 
426 void AliVZEROQADataMakerRec::InitDigits()
427 {
428   // create Digits histograms in Digits subdir
429   const Bool_t expert   = kTRUE ; 
430   const Bool_t image    = kTRUE ; 
431   
432   TH1I *fhDigTDC[64]; 
433   TH1I *fhDigADC[64]; 
434   
435   // create Digits histograms in Digits subdir
436   TH1I * h0 = new TH1I("hDigitMultiplicity", "Digits multiplicity distribution in VZERO;# of Digits;Entries", 100, 0, 99) ; 
437   h0->Sumw2() ;
438   Add2DigitsList(h0, 0, !expert, image) ;
439   
440   for (Int_t i=0; i<64; i++)
441     {
442       fhDigTDC[i] = new TH1I(Form("hDigitTDC%d", i),Form("Digit TDC in cell %d; TDC value;Entries",i),300,0.,149.);
443     
444       fhDigADC[i]= new TH1I(Form("hDigitADC%d",i),Form("Digit ADC in cell %d;ADC value;Entries",i),1024,0.,1023.);
445     
446       Add2DigitsList(fhDigTDC[i],i+1, !expert, image);
447       Add2DigitsList(fhDigADC[i],i+1+64, !expert, image);  
448     }  
449   //
450   ClonePerTrigClass(AliQAv1::kDIGITS); // this should be the last line
451 }
452
453 //____________________________________________________________________________
454 void AliVZEROQADataMakerRec::MakeDigits()
455 {
456   // makes data from Digits
457
458   FillDigitsData(0,fDigitsArray->GetEntriesFast()) ; 
459   TIter next(fDigitsArray) ; 
460   AliVZEROdigit *aVZERODigit ; 
461   while ( (aVZERODigit = dynamic_cast<AliVZEROdigit *>(next())) ) {
462     Int_t   aPMNumber  = aVZERODigit->PMNumber();         
463     FillDigitsData(aPMNumber +1,    aVZERODigit->Time()) ;    // in 100 of picoseconds
464     FillDigitsData(aPMNumber +1+64, aVZERODigit->ADC()) ;
465   }  
466 }
467
468
469 //____________________________________________________________________________
470 void AliVZEROQADataMakerRec::MakeDigits(TTree *digitTree)
471 {
472   // makes data from Digit Tree
473         
474   if ( fDigitsArray ) 
475     fDigitsArray->Clear() ; 
476   else 
477     fDigitsArray = new TClonesArray("AliVZEROdigit", 1000) ; 
478   
479   TBranch * branch = digitTree->GetBranch("VZERODigit") ;
480   if ( ! branch ) {
481     AliWarning("VZERO branch in Digit Tree not found") ; 
482   } else {
483     branch->SetAddress(&fDigitsArray) ;
484     branch->GetEntry(0) ; 
485     MakeDigits() ; 
486   }  
487   //  
488   IncEvCountCycleDigits();
489   IncEvCountTotalDigits();
490   //
491 }
492
493
494 //____________________________________________________________________________
495 void AliVZEROQADataMakerRec::MakeESDs(AliESDEvent * esd)
496 {
497   // Creates QA data from ESDs
498   
499   UInt_t eventType = esd->GetEventType();
500
501   switch (eventType){
502   case PHYSICS_EVENT:
503     AliESDVZERO *esdVZERO=esd->GetVZEROData();
504    
505     if (!esdVZERO) break;
506                   
507     FillESDsData(kCellMultiV0A,esdVZERO->GetNbPMV0A());
508     FillESDsData(kCellMultiV0C,esdVZERO->GetNbPMV0C());  
509     FillESDsData(kMIPMultiV0A,esdVZERO->GetMTotV0A());
510     FillESDsData(kMIPMultiV0C,esdVZERO->GetMTotV0C());  
511         
512     for(Int_t i=0;i<64;i++) {
513       FillESDsData(kMIPMultiChannel,(Float_t) i,(Float_t) esdVZERO->GetMultiplicity(i));
514       FillESDsData(kChargeChannel,(Float_t) i,(Float_t) esdVZERO->GetAdc(i));
515       if (i < 32) {
516         if(esdVZERO->BBTriggerV0C(i)) FillESDsData(kBBFlag,(Float_t) i);
517         if(esdVZERO->BGTriggerV0C(i)) FillESDsData(kBGFlag,(Float_t) i);
518       }
519       else {
520         if(esdVZERO->BBTriggerV0A(i-32)) FillESDsData(kBBFlag,(Float_t) i);  
521         if(esdVZERO->BGTriggerV0A(i-32)) FillESDsData(kBGFlag,(Float_t) i);
522       }                 
523       Float_t time = (Float_t) esdVZERO->GetTime(i); //Convert in ns:  1 TDC channel = 100ps 
524       FillESDsData(kTimeChannel,(Float_t) i,time);
525     }
526                                 
527     Float_t timeV0A = esdVZERO->GetV0ATime();
528     Float_t timeV0C = esdVZERO->GetV0CTime();
529     Float_t diffTime;
530
531     if(timeV0A<-1024.+1.e-6 || timeV0C<-1024.+1.e-6) diffTime = -1024.;
532     else diffTime = timeV0A - timeV0C;
533
534     FillESDsData(kESDV0ATime,timeV0A);
535     FillESDsData(kESDV0CTime,timeV0C);
536     FillESDsData(kESDDiffTime,diffTime);
537                 
538     break;
539   }  
540   //
541   IncEvCountCycleESDs();
542   IncEvCountTotalESDs();  
543   //
544 }
545
546 //____________________________________________________________________________
547 void AliVZEROQADataMakerRec::MakeRaws(AliRawReader* rawReader)
548 {
549   // Fills histograms with Raws, computes average ADC values dynamically (pedestal subtracted)
550                   
551                                           
552   // Check id histograms already created for this Event Specie
553   if ( ! GetRawsData(kPedestalInt0) )
554     InitRaws() ;
555
556   rawReader->Reset() ; 
557   AliVZERORawStream* rawStream  = new AliVZERORawStream(rawReader); 
558   if(!(rawStream->Next())) return;  
559  
560   eventTypeType eventType = rawReader->GetType();
561
562   Int_t    mulV0A = 0 ; 
563   Int_t    mulV0C = 0 ; 
564   Double_t timeV0A =0., timeV0C = 0.;
565   Double_t weightV0A =0., weightV0C = 0.;
566   UInt_t   itimeV0A=0, itimeV0C=0;
567   Double_t chargeV0A=0., chargeV0C=0.;
568   Double_t mipV0A=0., mipV0C=0.;
569
570   Double_t diffTime=-100000.;
571
572   
573   switch (eventType){
574   case PHYSICS_EVENT:
575   
576     //    fNTotEvents++; // Use framework counters instead
577
578     Int_t  iFlag=0;
579     Int_t  pedestal;
580     Int_t  integrator;
581     Bool_t flagBB[64];   
582     Bool_t flagBG[64];   
583     Int_t  mbCharge;
584     Float_t charge;
585     Int_t  offlineCh;
586     Float_t adc[64], time[64], width[64], timeCorr[64]; 
587
588     for(Int_t iChannel=0; iChannel<64; iChannel++) { // BEGIN : Loop over channels
589                    
590       offlineCh = rawStream->GetOfflineChannel(iChannel);
591                    
592       // Fill Pedestal histograms
593            
594       for(Int_t j=15; j<21; j++) {
595         if((rawStream->GetBGFlag(iChannel,j) || rawStream->GetBBFlag(iChannel,j))) iFlag++;
596       }
597
598       if(iFlag == 0){ //No Flag found
599         for(Int_t j=15; j<21; j++){
600           pedestal= (Int_t) rawStream->GetPedestal(iChannel, j);
601           integrator = rawStream->GetIntegratorFlag(iChannel, j);
602
603           FillRawsData((integrator == 0 ? kPedestalInt0 : kPedestalInt1),offlineCh,pedestal);
604         }
605       }
606
607       // Fill Charge EoI histograms
608            
609       adc[offlineCh]    = 0.0;
610
611       // Search for the maximum charge in the train of 21 LHC clocks 
612       // regardless of the integrator which has been operated:
613       Float_t maxadc = 0;
614       Int_t imax = -1;
615       Float_t adcPedSub[21];
616       for(Int_t iClock=0; iClock<21; iClock++){
617         Bool_t iIntegrator = rawStream->GetIntegratorFlag(iChannel,iClock);
618         Int_t k = offlineCh+64*iIntegrator;
619
620         //printf(Form("clock = %d adc = %f ped %f\n",iClock,rawStream->GetPedestal(iChannel,iClock),fPedestal[k]));
621
622         adcPedSub[iClock] = rawStream->GetPedestal(iChannel,iClock) - fCalibData->GetPedestal(k);
623         //                              if(adcPedSub[iClock] <= GetRecoParam()->GetNSigmaPed()*fCalibData->GetSigma(k)) {
624         if(adcPedSub[iClock] <= 2.*fCalibData->GetSigma(k)) {
625           adcPedSub[iClock] = 0;
626           continue;
627         }
628         //                              if(iClock < GetRecoParam()->GetStartClock() || iClock > GetRecoParam()->GetEndClock()) continue;
629         if(iClock < 8 || iClock > 12) continue;
630         if(adcPedSub[iClock] > maxadc) {
631           maxadc = adcPedSub[iClock];
632           imax   = iClock;
633         }
634       }
635       //printf(Form("Channel %d (online), %d (offline)\n",iChannel,j)); 
636       if (imax != -1) {
637         //                              Int_t start = imax - GetRecoParam()->GetNPreClocks();
638         Int_t start = imax - 2;
639         if (start < 0) start = 0;
640         //                              Int_t end = imax + GetRecoParam()->GetNPostClocks();
641         Int_t end = imax + 1;
642         if (end > 20) end = 20;
643         for(Int_t iClock = start; iClock <= end; iClock++) {
644           adc[offlineCh] += adcPedSub[iClock];
645         }
646       }
647         
648                 
649       Int_t iClock  = imax;
650       charge = rawStream->GetPedestal(iChannel,iClock); // Charge at the maximum 
651
652       integrator    = rawStream->GetIntegratorFlag(iChannel,iClock);
653       flagBB[offlineCh]  = rawStream->GetBBFlag(iChannel, iClock);
654       flagBG[offlineCh]  = rawStream->GetBGFlag(iChannel,iClock );
655       Int_t board = AliVZEROCalibData::GetBoardNumber(offlineCh);
656       time[offlineCh] = rawStream->GetTime(iChannel)*fCalibData->GetTimeResolution(board);
657       width[offlineCh] = rawStream->GetWidth(iChannel)*fCalibData->GetWidthResolution(board);
658
659       if (time[offlineCh] >= 1e-6) FillRawsData(kChargeEoI,offlineCh,adc[offlineCh]);
660
661       FillRawsData((integrator == 0 ? kChargeEoIInt0 : kChargeEoIInt1),offlineCh,charge);
662       if(flagBB[offlineCh]) FillRawsData((integrator == 0 ? kChargeEoIBBInt0 : kChargeEoIBBInt1),offlineCh,charge);
663       if(flagBG[offlineCh]) FillRawsData((integrator == 0 ? kChargeEoIBGInt0 : kChargeEoIBGInt1),offlineCh,charge);
664
665       Float_t sigma = fCalibData->GetSigma(offlineCh+64*integrator);
666
667                    
668       // Calculation of the number of MIP
669       Double_t mipEoI = adc[offlineCh] * fCalibData->GetMIPperADC(offlineCh);
670
671       if((adc[offlineCh] > 2.*sigma) && !(time[offlineCh] <1.e-6)){ 
672         FillRawsData(kRawMIPChannel,offlineCh,mipEoI);
673         if(offlineCh<32) {
674           mulV0C++;
675           chargeV0C += adc[offlineCh];
676           mipV0C += mipEoI;
677         } else {
678           mulV0A++;
679           chargeV0A += adc[offlineCh];
680           mipV0A += mipEoI;
681         }
682       }
683
684       // Fill Charge Minimum Bias Histograms
685                    
686       int idx;
687       for(Int_t iBunch=0; iBunch<10; iBunch++){
688         integrator = rawStream->GetIntMBFlag(iChannel, iBunch);
689         bool bbFlag     = rawStream->GetBBMBFlag(iChannel, iBunch);
690         bool bgFlag     = rawStream->GetBGMBFlag(iChannel, iBunch);
691         mbCharge   = rawStream->GetChargeMB(iChannel, iBunch);
692
693         if(integrator==0){
694           if(bbFlag==0){
695             if(bgFlag==0) idx = kChargeMBBB0BG0Int0;
696             else idx = kChargeMBBB0BG1Int0;
697           } else {
698             if(bgFlag==0) idx = kChargeMBBB1BG0Int0;
699             else idx = kChargeMBBB1BG1Int0;
700           }
701         } else {
702           if(bbFlag==0){
703             if(bgFlag==0) idx = kChargeMBBB0BG0Int1;
704             else idx = kChargeMBBB0BG1Int1;
705           } else {
706             if(bgFlag==0) idx = kChargeMBBB1BG0Int1;
707             else idx = kChargeMBBB1BG1Int1;
708           }
709         }
710         FillRawsData(idx,offlineCh,mbCharge);
711       }   
712
713       // Fill HPTDC Time Histograms
714       timeCorr[offlineCh] = CorrectLeadingTime(offlineCh,time[offlineCh],adc[offlineCh]);
715
716       const Float_t p1 = 2.50; // photostatistics term in the time resolution
717       const Float_t p2 = 3.00; // slewing related term in the time resolution
718       if(timeCorr[offlineCh]>-1024 + 1.e-6){
719         Float_t nphe = adc[offlineCh]*kChargePerADC/(fCalibData->GetGain(offlineCh)*TMath::Qe());
720         Float_t timeErr = 0;
721         if (nphe>1.e-6) timeErr = TMath::Sqrt(kIntTimeRes*kIntTimeRes+
722                                               p1*p1/nphe+
723                                               p2*p2*(fTimeSlewing->GetParameter(0)*fTimeSlewing->GetParameter(1))*(fTimeSlewing->GetParameter(0)*fTimeSlewing->GetParameter(1))*
724                                               TMath::Power(adc[offlineCh]/fCalibData->GetCalibDiscriThr(offlineCh,kTRUE),2.*(fTimeSlewing->GetParameter(1)-1.))/
725                                               (fCalibData->GetCalibDiscriThr(offlineCh,kTRUE)*fCalibData->GetCalibDiscriThr(offlineCh,kTRUE)));
726
727         if (timeErr>1.e-6) {
728           if (offlineCh<32) {
729             itimeV0C++;
730             timeV0C += timeCorr[offlineCh]/(timeErr*timeErr);
731             weightV0C += 1./(timeErr*timeErr);
732           }else{
733             itimeV0A++;
734             timeV0A += timeCorr[offlineCh]/(timeErr*timeErr);
735             weightV0A += 1./(timeErr*timeErr);
736           }
737         }
738       }
739       FillRawsData(kHPTDCTime,offlineCh,timeCorr[offlineCh]);
740       FillRawsData(kWidth,offlineCh,width[offlineCh]);
741       if(flagBB[offlineCh]) {
742         FillRawsData(kHPTDCTimeBB,offlineCh,timeCorr[offlineCh]);
743         FillRawsData(kWidthBB,offlineCh,width[offlineCh]);
744       }
745       if(flagBG[offlineCh]) {
746         FillRawsData(kHPTDCTimeBG,offlineCh,timeCorr[offlineCh]);
747         FillRawsData(kWidthBG,offlineCh,width[offlineCh]);
748       }
749
750       // Fill Flag and Charge Versus LHC-Clock histograms
751            
752       for(Int_t iEvent=0; iEvent<21; iEvent++){
753         charge = rawStream->GetPedestal(iChannel,iEvent);
754         integrator = rawStream->GetIntegratorFlag(iChannel,iEvent);
755         bool bbFlag       = rawStream->GetBBFlag(iChannel, iEvent);
756         bool bgFlag       = rawStream->GetBGFlag(iChannel,iEvent );
757
758         FillRawsData((integrator == 0 ? kChargeVsClockInt0 : kChargeVsClockInt1 ), offlineCh,(float)iEvent-10,(float)charge);
759         FillRawsData(kBBFlagVsClock, offlineCh,(float)iEvent-10,(float)bbFlag);
760         FillRawsData(kBGFlagVsClock, offlineCh,(float)iEvent-10,(float)bgFlag);
761         if(iEvent==10) FillRawsData(kBBFlagsPerChannel, offlineCh,(float)bbFlag);
762       }
763
764     }// END of Loop over channels
765
766     if(weightV0A>1.e-6) timeV0A /= weightV0A; 
767     else timeV0A = -1024.;
768     if(weightV0C>1.e-6) timeV0C /= weightV0C;
769     else timeV0C = -1024.;
770     if(timeV0A<-1024.+1.e-6 || timeV0C<-1024.+1.e-6) diffTime = -1024.;
771     else diffTime = timeV0A - timeV0C;
772     
773     Bool_t v0ABB = kFALSE;
774     Bool_t v0CBB = kFALSE;
775     Bool_t v0ABG = kFALSE;
776     Bool_t v0CBG = kFALSE;
777     Bool_t v0AFake = kFALSE;
778     Bool_t v0CFake = kFALSE;
779     Bool_t v0AEmpty = kFALSE;
780     Bool_t v0CEmpty = kFALSE;
781     Int_t v0ATrigger=0;
782     Int_t v0CTrigger=0; 
783     
784     // Change default BB and BG windows according to the Trigger Count Offset setting with respect to the default one which is 3247.
785     Float_t winOffset = (fCalibData->GetTriggerCountOffset(0) - 3247)*25.;
786     
787     if((timeV0A>kMinBBA-winOffset) && (timeV0A<kMaxBBA-winOffset)) {
788       v0ABB = kTRUE;
789       v0ATrigger=2;
790     } else if((timeV0A>kMinBGA-winOffset) && (timeV0A<kMaxBGA-winOffset)) {
791       v0ABG = kTRUE;
792       v0ATrigger=3;
793     } else if(timeV0A>-1024.+1.e-6) {
794       v0AFake = kTRUE;
795       v0ATrigger=1;
796     } else {
797       v0AEmpty = kTRUE;
798       v0ATrigger=0;
799     }
800     
801     if((timeV0C>kMinBBC-winOffset) && (timeV0C<kMaxBBC-winOffset)) {
802       v0CBB = kTRUE;
803       v0CTrigger=2;
804     } else if((timeV0C>kMinBGC-winOffset) && (timeV0C<kMaxBGC-winOffset)) {
805       v0CBG = kTRUE;
806       v0CTrigger=3;
807     } else if(timeV0C>-1024.+1.e-6) {
808       v0CFake = kTRUE;
809       v0CTrigger=1;
810     } else {
811       v0CEmpty = kTRUE;
812       v0CTrigger=0;
813     }
814
815     // Fill Trigger output histogram
816     if(v0ABB && v0CBB) FillRawsData(kTriggers,0);
817     if((v0ABB || v0CBB) && !(v0ABG || v0CBG)) FillRawsData(kTriggers,1);
818     if(v0ABG && v0CBB) FillRawsData(kTriggers,2);
819     if(v0ABB && v0CBG) FillRawsData(kTriggers,3);
820                 
821     FillRawsData(kTriggers2,v0ATrigger,v0CTrigger);
822
823     FillRawsData(kV0ATime,timeV0A);
824     FillRawsData(kV0CTime,timeV0C);
825     FillRawsData(kDiffTime,diffTime);
826     FillRawsData(kTimeV0AV0C,timeV0A,timeV0C);
827
828     FillRawsData(kMultiV0A,mulV0A);
829     FillRawsData(kMultiV0C,mulV0C);
830
831     FillRawsData(kChargeV0A,chargeV0A);
832     FillRawsData(kChargeV0C,chargeV0C);
833     FillRawsData(kChargeV0,chargeV0A + chargeV0C);
834
835     FillRawsData(kRawMIPV0A,mipV0A);
836     FillRawsData(kRawMIPV0C,mipV0C);
837     FillRawsData(kRawMIPV0,mipV0A + mipV0C);
838     break;
839             
840   } // END of SWITCH : EVENT TYPE 
841         
842   //  fEvent++;  // RS: Use framework counters instead
843   TParameter<double> * p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kMultiV0A)->GetName()))) ; 
844   if (p) p->SetVal((double)mulV0A) ; 
845
846   p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kMultiV0C)->GetName()))) ; 
847   if (p) p->SetVal((double)mulV0C) ;                     
848
849   p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kChargeV0A)->GetName()))) ; 
850   if (p) p->SetVal((double)chargeV0A) ; 
851
852   p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kChargeV0C)->GetName()))) ; 
853   if (p) p->SetVal((double)chargeV0C) ;                     
854
855   p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kChargeV0)->GetName()))) ; 
856   if (p) p->SetVal((double)(chargeV0A + chargeV0C)) ;                     
857         
858   p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kRawMIPV0A)->GetName()))) ; 
859   if (p) p->SetVal((double)mipV0A) ; 
860         
861   p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kRawMIPV0C)->GetName()))) ; 
862   if (p) p->SetVal((double)mipV0C) ;                     
863         
864   p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kRawMIPV0)->GetName()))) ; 
865   if (p) p->SetVal((double)(mipV0A + mipV0C)) ;                     
866         
867   p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kV0ATime)->GetName()))) ; 
868   if (p) p->SetVal((double)timeV0A) ; 
869         
870   p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kV0CTime)->GetName()))) ; 
871   if (p) p->SetVal((double)timeV0C) ;                     
872         
873   p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kDiffTime)->GetName()))) ; 
874   if (p) p->SetVal((double)diffTime) ;                     
875         
876   delete rawStream; rawStream = 0x0;      
877   //
878   IncEvCountCycleRaws();
879   IncEvCountTotalRaws();
880   //
881 }
882
883 //____________________________________________________________________________ 
884 void AliVZEROQADataMakerRec::StartOfDetectorCycle()
885 {
886   // Detector specific actions at start of cycle
887   
888   // Reset of the histogram used - to have the trend versus time -
889  
890   fCalibData = GetCalibData();
891  
892   AliCDBEntry *entry = AliCDBManager::Instance()->Get("GRP/CTP/CTPtiming");
893   if (!entry) AliFatal("CTP timing parameters are not found in OCDB !");
894   AliCTPTimeParams *ctpParams = (AliCTPTimeParams*)entry->GetObject();
895   Float_t l1Delay = (Float_t)ctpParams->GetDelayL1L0()*25.0;
896
897   AliCDBEntry *entry1 = AliCDBManager::Instance()->Get("GRP/CTP/TimeAlign");
898   if (!entry1) AliFatal("CTP time-alignment is not found in OCDB !");
899   AliCTPTimeParams *ctpTimeAlign = (AliCTPTimeParams*)entry1->GetObject();
900   l1Delay += ((Float_t)ctpTimeAlign->GetDelayL1L0()*25.0);
901
902   AliCDBEntry *entry2 = AliCDBManager::Instance()->Get("VZERO/Calib/TimeDelays");
903   if (!entry2) AliFatal("VZERO time delays are not found in OCDB !");
904   TH1F *delays = (TH1F*)entry2->GetObject();
905
906   AliCDBEntry *entry3 = AliCDBManager::Instance()->Get("VZERO/Calib/TimeSlewing");
907   if (!entry3) AliFatal("VZERO time slewing function is not found in OCDB !");
908   fTimeSlewing = (TF1*)entry3->GetObject();
909
910   for(Int_t i = 0 ; i < 64; ++i) {
911     //Int_t board = AliVZEROCalibData::GetBoardNumber(i);
912     fTimeOffset[i] = (
913                       //        ((Float_t)fCalibData->GetTriggerCountOffset(board) -
914                       //        (Float_t)fCalibData->GetRollOver(board))*25.0 +
915                       //     fCalibData->GetTimeOffset(i) -
916                       //     l1Delay+
917                       delays->GetBinContent(i+1)//+
918                       //      kV0Offset
919                       );
920     //                AliInfo(Form(" fTimeOffset[%d] = %f  kV0offset %f",i,fTimeOffset[i],kV0Offset));
921   }
922
923  
924  
925   
926         
927   TTimeStamp currentTime;
928   fCycleStartTime = currentTime.GetSec();
929  
930   //  fNTotEvents = 0;
931 }
932
933
934 //-------------------------------------------------------------------------------------------------
935 Float_t AliVZEROQADataMakerRec::CorrectLeadingTime(Int_t i, Float_t time, Float_t adc) const
936 {
937   // Correct the leading time
938   // for slewing effect and
939   // misalignment of the channels
940   if (time < 1e-6) return -1024;
941
942   // Channel alignment and general offset subtraction
943   //  if (i < 32) time -= kV0CDelayCables;
944   //  time -= fTimeOffset[i];
945   //AliInfo(Form("time-offset %f", time));
946
947   // In case of pathological signals
948   if (adc < 1e-6) return time;
949
950   // Slewing correction
951   Float_t thr = fCalibData->GetCalibDiscriThr(i,kTRUE);
952   //AliInfo(Form("adc %f thr %f dtime %f ", adc,thr,fTimeSlewing->Eval(adc/thr)));
953   time -= fTimeSlewing->Eval(adc/thr);
954
955   return time;
956 }
957