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