]> git.uio.no Git - u/mrichter/AliRoot.git/blob - VZERO/AliVZEROQADataMakerRec.cxx
Adding optional check on the consistency of the ALTRO payload. By default it is on...
[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 <TH1F.h> 
26 #include <TH1I.h> 
27 #include <TH2I.h> 
28 #include <TH2D.h> 
29 #include <TGraph.h> 
30 #include <TParameter.h>
31
32 // --- Standard library ---
33
34 // --- AliRoot header files ---
35 #include "AliESDEvent.h"
36 #include "AliLog.h"
37 #include "AliCDBManager.h"
38 #include "AliCDBStorage.h"
39 #include "AliCDBEntry.h"
40 #include "AliVZEROQADataMakerRec.h"
41 #include "AliQAChecker.h"
42 #include "AliRawReader.h"
43 #include "AliVZERORawStream.h"
44 #include "AliVZEROdigit.h"
45 #include "AliVZEROReconstructor.h"
46 #include "event.h"
47
48
49 ClassImp(AliVZEROQADataMakerRec)
50            
51 //____________________________________________________________________________ 
52   AliVZEROQADataMakerRec::AliVZEROQADataMakerRec() : 
53         AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kVZERO), "VZERO Quality Assurance Data Maker"),
54         fCalibData(0x0),
55     fEvent(0)
56     
57 {
58    // Constructor
59    
60       AliDebug(AliQAv1::GetQADebugLevel(), "Construct VZERO QA Object");
61   
62    for(Int_t i=0; i<64; i++){  
63        fEven[i] = 0;   
64        fOdd[i]  = 0;  }
65   
66    for(Int_t i=0; i<128; i++){  
67        fADCmean[i] = 0.0;   }   
68 }
69
70 //____________________________________________________________________________ 
71   AliVZEROQADataMakerRec::AliVZEROQADataMakerRec(const AliVZEROQADataMakerRec& qadm) :
72   AliQADataMakerRec(),
73         fCalibData(0x0),
74     fEvent(0)
75   
76 {
77   // Copy constructor 
78   
79    SetName((const char*)qadm.GetName()) ; 
80    SetTitle((const char*)qadm.GetTitle()); 
81 }
82
83 //__________________________________________________________________
84 AliVZEROQADataMakerRec& AliVZEROQADataMakerRec::operator = (const AliVZEROQADataMakerRec& qadm )
85 {
86   // Equal operator
87   
88   this->~AliVZEROQADataMakerRec();
89   new(this) AliVZEROQADataMakerRec(qadm);
90   return *this;
91 }
92
93 //____________________________________________________________________________
94 AliVZEROCalibData* AliVZEROQADataMakerRec::GetCalibData() const
95
96 {
97    AliCDBManager *man = AliCDBManager::Instance();
98
99    AliCDBEntry *entry=0;
100
101    entry = man->Get("VZERO/Calib/Data",fRun);
102    if(!entry){
103         AliWarning("Load of calibration data from default storage failed!");
104         AliWarning("Calibration data will be loaded from local storage ($ALICE_ROOT)");
105         
106         man->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
107         entry = man->Get("VZERO/Calib/Data",fRun);
108    }
109    // Retrieval of data in directory VZERO/Calib/Data:
110
111    AliVZEROCalibData *calibdata = 0;
112
113    if (entry) calibdata = (AliVZEROCalibData*) entry->GetObject();
114    if (!calibdata)  AliFatal("No calibration data from calibration database !");
115
116    return calibdata;
117 }
118
119
120  
121 //____________________________________________________________________________ 
122 void AliVZEROQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
123 {
124   // Detector specific actions at end of cycle
125   // Does the QA checking
126   
127   AliQAChecker::Instance()->Run(AliQAv1::kVZERO, task, list) ;
128
129   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
130     if (! IsValidEventSpecie(specie, list)) 
131       continue ;
132     SetEventSpecie(AliRecoParam::ConvertIndex(specie)) ; 
133     if(task == AliQAv1::kRAWS){
134           int nMaxBin = GetRawsData(kPedestalTimeInt0)->GetNbinsY();
135       if(fCurrentCycle%nMaxBin==0) {
136         GetRawsData(kPedestalTimeInt0)->Reset();
137         GetRawsData(kPedestalTimeInt1)->Reset();
138         GetRawsData(kChargeEoITimeInt0)->Reset();
139         GetRawsData(kChargeEoITimeInt1)->Reset();
140       }
141       TH1D* hProj;
142       char name[50];
143       for(Int_t iChannel=0; iChannel<64; iChannel++) {
144         for(Int_t integrator=0;integrator<2;integrator++){
145           sprintf(name,"Ped_%d_%d",iChannel,integrator);
146           hProj = ((TH2I*)GetRawsData((integrator == 0 ? kPedestalCycleInt0 : kPedestalCycleInt1)))->ProjectionY(name,iChannel+1,iChannel+1);
147           ((TH2D*)GetRawsData((integrator == 0 ? kPedestalTimeInt0 : kPedestalTimeInt1)))->Fill((double)iChannel,(double)(fCurrentCycle%nMaxBin),(double)hProj->GetMean());
148           delete hProj;
149
150           sprintf(name,"Charge_%d_%d",iChannel,integrator);
151           hProj = ((TH2I*)GetRawsData((integrator == 0 ? kChargeEoICycleInt0 : kChargeEoICycleInt1)))->ProjectionY(name,iChannel+1,iChannel+1);
152           ((TH2D*)GetRawsData((integrator == 0 ? kChargeEoITimeInt0 : kChargeEoITimeInt1)))->Fill((double)iChannel,(double)(fCurrentCycle%nMaxBin),hProj->GetMean());
153           delete hProj;
154         }
155       }
156     } else if (task == AliQAv1::kESDS) {
157     }
158   }
159 }
160
161 //____________________________________________________________________________ 
162 void AliVZEROQADataMakerRec::InitESDs()
163 {
164   // Creates histograms to control ESDs
165   
166   const Bool_t expert   = kTRUE ; 
167   const Bool_t image    = kTRUE ; 
168         
169   TH2D * h2d;
170   TH1I * h1i;
171   TH1D * h1d;
172                 
173   h1i = new TH1I("H1I_Cell_Multiplicity_V0A", "Cell Multiplicity in V0A;Multiplicity (Nb of Cell);Counts", 35, 0, 35) ;  
174   Add2ESDsList(h1i, kCellMultiV0A, !expert, image)  ;  
175                                                                                                         
176   h1i = new TH1I("H1I_Cell_Multiplicity_V0C", "Cell Multiplicity in V0;Multiplicity (Nb of Cell);Counts", 35, 0, 35) ;  
177   Add2ESDsList(h1i, kCellMultiV0C, !expert, image)  ;  
178    
179   h1d = new TH1D("H1D_MIP_Multiplicity_V0A", "MIP Multiplicity in V0A;Multiplicity (Nb of MIP);Counts", 1000, 0, 1000) ;  
180   Add2ESDsList(h1d, kMIPMultiV0A, !expert, image)  ;  
181   
182   h1d = new TH1D("H1D_MIP_Multiplicity_V0C", "MIP Multiplicity in V0C;Multiplicity (Nb of MIP);Counts", 1000, 0, 1000) ;  
183   Add2ESDsList(h1d, kMIPMultiV0C, !expert, image)  ;  
184
185   h2d = new TH2D("H2D_MIP_Multiplicity_Channel", "MIP Multiplicity per Channel;Channel;Multiplicity (Nb of MIP)",64, 0, 64, 100, 0, 100) ;  
186   Add2ESDsList(h2d, kMIPMultiChannel, !expert, image)  ;  
187   
188   h1d = new TH1D("H1D_BBFlag_Counters", "BB Flag Counters;Channel;Counts",64, 0, 64) ;  
189   Add2ESDsList(h1d, kBBFlag, !expert, image)  ;  
190   
191   h1d = new TH1D("H1D_BGFlag_Counters", "BG Flag Counters;Channel;Counts",64, 0, 64) ;  
192   Add2ESDsList(h1d, kBGFlag, !expert, image)  ;  
193   
194   h2d = new TH2D("H2D_Charge_Channel", "ADC Charge per channel;Channel;Charge (ADC counts)",64, 0, 64, 1024, 0, 1024) ;  
195   Add2ESDsList(h2d, kChargeChannel, !expert, image)  ;  
196   
197   h2d = new TH2D("H2D_Time_Channel", "Time per channel;Channel;Time (ns)",64, 0, 64, 820, 0, 410) ;  
198   Add2ESDsList(h2d, kTimeChannel, !expert, image)  ;  
199   
200   h1d = new TH1D("H1D_V0A_Time", "Mean V0A Time;Time (ns);Counts",2048, 0., 409.6);
201   Add2ESDsList(h1d,kESDV0ATime, !expert, image); 
202   
203   h1d = new TH1D("H1D_V0C_Time", "Mean V0C Time;Time (ns);Counts",2048, 0., 409.6);
204   Add2ESDsList(h1d,kESDV0CTime, !expert, image); 
205   
206   h1d = new TH1D("H1D_Diff_Time", "Diff Time V0A - V0C;Diff Time V0A - V0C (ns);Counts",2*2048, -409.6, 409.6);
207   Add2ESDsList(h1d,kESDDiffTime, !expert, image); 
208         
209 }
210
211 //____________________________________________________________________________ 
212  void AliVZEROQADataMakerRec::InitRaws()
213  {
214    // Creates RAW histograms in Raws subdir
215    
216    const Bool_t expert   = kTRUE ; 
217    const Bool_t saveCorr = kTRUE ; 
218    const Bool_t image    = kTRUE ; 
219
220   char name[50] , title[100];
221   const Int_t kNintegrator  =    2;
222  
223   const Int_t kNTdcTimeBins  = 2048;
224   const Int_t kTdcTimeMin    =    0;
225   const Int_t kTdcTimeMax    = 4096;
226   const Int_t kNTdcWidthBins =  128;
227   const Int_t kTdcWidthMin   =    0;
228   const Int_t kTdcWidthMax   =  128;
229   const Int_t kNChargeBins   = 1024;
230   const Int_t kChargeMin     =    0;
231   const Int_t kChargeMax     = 1024;
232   const Int_t kNChannelBins  =   64;
233   const Int_t kChannelMin    =    0;
234   const Int_t kChannelMax    =   64;
235   const Int_t kNPedestalBins =  200;
236   const Int_t kPedestalMin   =    0;
237   const Int_t kPedestalMax   =  200;
238   const Int_t kTimeMin       =   0;
239   const Int_t kTimeMax       = 100;
240   const Int_t kNMIPBins      = 200;
241   const Int_t kMIPMin        =   0;
242   const Int_t kMIPMax        = 200;
243
244   TH2I * h2i;
245   TH2D * h2d;
246   TH1I * h1i;
247   TH1D * h1d;
248
249   int iHisto =0;
250  
251    // Creation of Cell Multiplicity Histograms
252   h1i = new TH1I("H1I_Multiplicity_V0A", "Cell Multiplicity in V0A;# of Cells;Entries", 35, 0, 35) ;  
253   Add2RawsList(h1i,kMultiV0A, !expert, image, saveCorr);   iHisto++;
254   h1i = new TH1I("H1I_Multiplicity_V0C", "Cell Multiplicity in V0C;# of Cells;Entries", 35, 0, 35) ;  
255   Add2RawsList(h1i,kMultiV0C, !expert, image, saveCorr);   iHisto++;
256  
257   // Creation of Total Charge Histograms
258   h1d = new TH1D("H1D_Charge_V0A", "Total Charge in V0A;Charge [ADC counts];Counts", 2048, 0, 32768) ;  
259   Add2RawsList(h1d,kChargeV0A, !expert, image, saveCorr);   iHisto++;
260   h1d = new TH1D("H1D_Charge_V0C", "Total Charge in V0C;Charge [ADC counts];Counts", 2048, 0, 32768) ;  
261   Add2RawsList(h1d,kChargeV0C, !expert, image, saveCorr);   iHisto++;
262   h1d = new TH1D("H1D_Charge_V0", "Total Charge in V0;Charge [ADC counts];Counts", 2048, 0, 65536) ;  
263   Add2RawsList(h1d,kChargeV0, !expert, image, saveCorr);   iHisto++;
264   
265   // Creation of MIP Histograms
266   h1d = new TH1D("H1D_MIP_V0A", "Total MIP in V0A;Charge [MIP];Counts", 2*kNMIPBins,kMIPMin ,32*kMIPMax) ;  
267   Add2RawsList(h1d,kRawMIPV0A, !expert, image, saveCorr);   iHisto++;
268   h1d = new TH1D("H1D_MIP_V0C", "Total MIP in V0C;Charge [MIP];Counts", 2*kNMIPBins,kMIPMin ,32*kMIPMax) ;  
269   Add2RawsList(h1d,kRawMIPV0C, !expert, image, saveCorr);   iHisto++;
270   h1d = new TH1D("H1D_MIP_V0", "Total MIP in V0;Charge [MIP];Counts", 2*kNMIPBins,kMIPMin ,32*kMIPMax) ;  
271   Add2RawsList(h1d,kRawMIPV0, !expert, image, saveCorr);   iHisto++;
272   h2d = new TH2D("H2D_MIP_Channel", "Nb of MIP per channel;Channel;# of Mips", kNChannelBins, kChannelMin, kChannelMax,kNMIPBins,kMIPMin ,kMIPMax) ;  
273   Add2RawsList(h2d,kRawMIPChannel, expert, !image, !saveCorr);   iHisto++;
274   
275  
276  for(Int_t iInt=0;iInt<kNintegrator;iInt++){
277     // Creation of Pedestal histograms 
278     sprintf(name,"H2I_Pedestal_Int%d",iInt);
279     sprintf(title,"Pedestal (Int%d);Pedestal [ADC counts];Counts",iInt);
280     h2i = new TH2I(name, title,kNChannelBins, kChannelMin, kChannelMax,kNPedestalBins,kPedestalMin ,kPedestalMax );
281     Add2RawsList(h2i,(iInt == 0 ? kPedestalInt0 : kPedestalInt1), expert, !image, !saveCorr); iHisto++;
282         
283     // Creation of temporary Pedestal histo used for the mean versus time histogram. This histogram will be reset at the end of each cycle
284     sprintf(name,"H2I_Pedestal_CycleInt%d",iInt);
285     sprintf(title,"One Cycle Pedestal (Int%d);Pedestal [ADC counts];Counts",iInt);
286     h2i = new TH2I(name, title,kNChannelBins, kChannelMin, kChannelMax,kNPedestalBins,kPedestalMin ,kPedestalMax );
287     Add2RawsList(h2i,(iInt == 0 ? kPedestalCycleInt0 : kPedestalCycleInt1), expert, !image, !saveCorr); iHisto++;
288                 
289     // Creation of Pedestal versus time graph.
290     sprintf(name,"H2D_Pedestal_Time_Int%d",iInt);
291     sprintf(title,"Pedestal Versus Time (Int%d);Time [ns];Pedestal [ADC counts]",iInt);
292     h2d = new TH2D(name, title,kNChannelBins, kChannelMin, kChannelMax,kTimeMax,kTimeMin ,kTimeMax );
293     Add2RawsList(h2d,(iInt == 0 ? kPedestalTimeInt0 : kPedestalTimeInt1), expert, !image, !saveCorr); iHisto++;
294
295    // Creation of Charge EoI histograms 
296     sprintf(name,"H2I_ChargeEoI_Int%d",iInt);
297     sprintf(title,"Charge EoI (Int%d);Charge [ADC counts];Counts",iInt);
298     h2i = new TH2I(name, title,kNChannelBins, kChannelMin, kChannelMax, kNChargeBins, kChargeMin, kChargeMax);
299     Add2RawsList(h2i,(iInt == 0 ? kChargeEoIInt0 : kChargeEoIInt1), !expert, image, !saveCorr); iHisto++;
300
301    // Creation of temporary Charge EoI histograms used for the mean versus time histogram. This histogram will be reset at the end of each cycle
302     sprintf(name,"H2I_ChargeEoI_CycleInt%d",iInt);
303     sprintf(title,"One Cycle Charge EoI (Int%d);Charge [ADC counts];Counts",iInt);
304     h2i = new TH2I(name, title,kNChannelBins, kChannelMin, kChannelMax, kNChargeBins, kChargeMin, kChargeMax);
305     Add2RawsList(h2i,(iInt == 0 ? kChargeEoICycleInt0 : kChargeEoICycleInt1), expert, !image, !saveCorr); iHisto++;
306                 
307     // Creation of Charge EoI versus time graphs
308     sprintf(name,"H2D_ChargeEoI_Time_Int%d",iInt);
309     sprintf(title,"Charge EoI Versus Time (Int%d);Time [ns];Charge [ADC counts]",iInt);
310     h2d = new TH2D(name, title,kNChannelBins, kChannelMin, kChannelMax,kTimeMax,kTimeMin ,kTimeMax );
311     Add2RawsList(h2d,(iInt == 0 ? kChargeEoITimeInt0 : kChargeEoITimeInt1), expert, !image, !saveCorr); iHisto++;
312     
313     sprintf(name,"H2I_ChargeEoI_BB_Int%d",iInt);
314     sprintf(title,"Charge EoI w/ BB Flag (Int%d);??;Charge [ADC counts]",iInt);
315     h2i = new TH2I(name, title,kNChannelBins, kChannelMin, kChannelMax, kNChargeBins, kChargeMin, kChargeMax);
316     Add2RawsList(h2i,(iInt == 0 ? kChargeEoIBBInt0 : kChargeEoIBBInt1), expert, !image, !saveCorr); iHisto++;
317     
318     sprintf(name,"H2I_ChargeEoI_BG_Int%d",iInt);
319     sprintf(title,"Charge EoI w/ BG Flag (Int%d);??;Charge [ADC counts]",iInt);
320     h2i = new TH2I(name, title,kNChannelBins, kChannelMin, kChannelMax, kNChargeBins, kChargeMin, kChargeMax);
321     Add2RawsList(h2i,(iInt == 0 ?  kChargeEoIBGInt0: kChargeEoIBGInt1), expert, !image, !saveCorr); iHisto++;
322
323     // Creation of Charge versus LHC Clock histograms 
324     sprintf(name,"H2D_ChargeVsClock_Int%d",iInt);
325     sprintf(title,"Charge Versus LHC-Clock (Int%d);Tine [ns];Charge [ADC counts]",iInt);
326     h2d = new TH2D(name, title,kNChannelBins, kChannelMin, kChannelMax,21, -10.5, 10.5 );
327     Add2RawsList(h2d,(iInt == 0 ? kChargeVsClockInt0 : kChargeVsClockInt1 ), expert, !image, !saveCorr); iHisto++;
328         
329     // Creation of Minimum Bias Charge histograms 
330     for(Int_t iBB=0;iBB<2;iBB++){
331                 for(Int_t iBG=0;iBG<2;iBG++){
332                         sprintf(name,"H2I_ChargeMB_BB%d_BG%d_Int%d",iBB,iBG,iInt);
333                         sprintf(title,"MB Charge (BB=%d, BG=%d, Int=%d);Charge [ADC counts];Counts",iBB,iBG,iInt);
334                         h2i = new TH2I(name, title,kNChannelBins, kChannelMin, kChannelMax,kNChargeBins, kChargeMin, kChargeMax);
335                         int idx;
336                         if(iInt==0){
337                                 if(iBB==0){
338                                         if(iBG==0) idx = kChargeMBBB0BG0Int0;
339                                         else idx = kChargeMBBB0BG1Int0;
340                                 } else {
341                                         if(iBG==0) idx = kChargeMBBB1BG0Int0;
342                                         else idx = kChargeMBBB1BG1Int0;
343                                 }
344                         } else {
345                                 if(iBB==0){
346                                         if(iBG==0) idx = kChargeMBBB0BG0Int1;
347                                         else idx = kChargeMBBB0BG1Int1;
348                                 } else {
349                                         if(iBG==0) idx = kChargeMBBB1BG0Int1;
350                                         else idx = kChargeMBBB1BG1Int1;
351                                 }
352                         }
353                         Add2RawsList(h2i,idx, expert, !image, !saveCorr); iHisto++;
354                 }
355     }
356         
357  }
358  
359      // Creation of Time histograms 
360         sprintf(name,"H2I_Width");
361         sprintf(title,"HPTDC Width;Width [ns];Counts");
362         h2i = new TH2I(name, title,kNChannelBins, kChannelMin, kChannelMax, kNTdcWidthBins, kTdcWidthMin, kTdcWidthMax);
363         Add2RawsList(h2i,kWidth, expert, !image, !saveCorr); iHisto++;
364
365         sprintf(name,"H2I_Width_BB");
366         sprintf(title,"HPTDC Width w/ BB Flag condition;??;Width [ns]");
367         h2i = new TH2I(name, title,kNChannelBins, kChannelMin, kChannelMax, kNTdcWidthBins, kTdcWidthMin, kTdcWidthMax);
368         Add2RawsList(h2i,kWidthBB, expert, !image, !saveCorr); iHisto++;
369
370         sprintf(name,"H2I_Width_BG");
371         sprintf(title,"HPTDC Width w/ BG Flag condition??;Width [ns]");
372         h2i = new TH2I(name, title,kNChannelBins, kChannelMin, kChannelMax, kNTdcWidthBins, kTdcWidthMin, kTdcWidthMax);
373         Add2RawsList(h2i,kWidthBG, expert, !image, !saveCorr); iHisto++;
374
375         sprintf(name,"H2I_HPTDCTime");
376         sprintf(title,"HPTDC Time;??;Width [ns]");
377         h2i = new TH2I(name, title,kNChannelBins, kChannelMin, kChannelMax, kNTdcTimeBins, kTdcTimeMin, kTdcTimeMax);
378         Add2RawsList(h2i,kHPTDCTime, !expert, image, !saveCorr); iHisto++;
379
380         sprintf(name,"H2I_HPTDCTime_BB");
381         sprintf(title,"HPTDC Time w/ BB Flag condition;??;Width [ns]");
382         h2i = new TH2I(name, title,kNChannelBins, kChannelMin, kChannelMax, kNTdcTimeBins, kTdcTimeMin, kTdcTimeMax);
383         Add2RawsList(h2i,kHPTDCTimeBB, expert, !image, !saveCorr); iHisto++;
384
385         sprintf(name,"H2I_HPTDCTime_BG");
386         sprintf(title,"HPTDC Time w/ BG Flag condition;??;Width [ns]");
387         h2i = new TH2I(name, title,kNChannelBins, kChannelMin, kChannelMax, kNTdcTimeBins, kTdcTimeMin, kTdcTimeMax);
388         Add2RawsList(h2i,kHPTDCTimeBG, expert, !image, !saveCorr); iHisto++;
389         
390         sprintf(name,"H1D_V0A_Time");
391         sprintf(title,"V0A Time;Time [ns];Counts");
392         h1d = new TH1D(name, title,kNTdcTimeBins, kTdcTimeMin/10, kTdcTimeMax/10);
393         Add2RawsList(h1d,kV0ATime, !expert, !image, saveCorr); iHisto++;
394         
395         sprintf(name,"H1D_V0C_Time");
396         sprintf(title,"V0C Time;Time [ns];Counts");
397         h1d = new TH1D(name, title,kNTdcTimeBins, kTdcTimeMin/10, kTdcTimeMax/10);
398         Add2RawsList(h1d,kV0CTime, !expert, !image, saveCorr); iHisto++;
399         
400         sprintf(name,"H1D_Diff_Time");
401         sprintf(title,"Diff V0A-V0C Time;Time [ns];Counts");
402         h1d = new TH1D(name, title,2*kNTdcTimeBins, -kTdcTimeMax/10, kTdcTimeMax/10);
403         Add2RawsList(h1d,kDiffTime, !expert, !image, saveCorr); iHisto++;
404         
405         // Creation of Flag versus LHC Clock histograms 
406         sprintf(name,"H2D_BBFlagVsClock");
407         sprintf(title,"BB-Flags Versus LHC-Clock;Time [ns];??");
408         h2d = new TH2D(name, title,kNChannelBins, kChannelMin, kChannelMax,21, -10.5, 10.5 );
409         Add2RawsList(h2d,kBBFlagVsClock, expert, !image, !saveCorr); iHisto++;
410         
411         sprintf(name,"H2D_BGFlagVsClock");
412         sprintf(title,"BG-Flags Versus LHC-Clock;Time [ns];??");
413         h2d = new TH2D(name, title,kNChannelBins, kChannelMin, kChannelMax,21, -10.5, 10.5 );
414         Add2RawsList(h2d,kBGFlagVsClock, expert, !image, !saveCorr); iHisto++;
415          
416         AliDebug(AliQAv1::GetQADebugLevel(), Form("%d Histograms has been added to the Raws List",iHisto));
417  }
418
419 //____________________________________________________________________________ 
420 void AliVZEROQADataMakerRec::InitDigits()
421 {
422   // create Digits histograms in Digits subdir
423   const Bool_t expert   = kTRUE ; 
424   const Bool_t image    = kTRUE ; 
425   
426   char TDCname[100];
427   char ADCname[100];
428   TH1I *fhDigTDC[64]; 
429   TH1I *fhDigADC[64]; 
430   char texte[100];
431   
432   // create Digits histograms in Digits subdir
433   TH1I * h0 = new TH1I("hDigitMultiplicity", "Digits multiplicity distribution in VZERO;# of Digits;Entries", 100, 0, 99) ; 
434   h0->Sumw2() ;
435   Add2DigitsList(h0, 0, !expert, image) ;
436   
437   for (Int_t i=0; i<64; i++)
438     {
439     sprintf(TDCname, "hDigitTDC%d", i);
440     sprintf(texte,"Digit TDC in cell %d; TDC value;Entries",i);    
441     fhDigTDC[i] = new TH1I(TDCname,texte,300,0.,149.);
442     
443     sprintf(ADCname,"hDigitADC%d",i);
444     sprintf(texte,"Digit ADC in cell %d;ADC value;Entries",i);
445     fhDigADC[i]= new TH1I(ADCname,texte,1024,0.,1023.);
446     
447     Add2DigitsList(fhDigTDC[i],i+1, !expert, image);
448     Add2DigitsList(fhDigADC[i],i+1+64, !expert, image);  
449     }  
450 }
451
452 //____________________________________________________________________________
453 void AliVZEROQADataMakerRec::MakeDigits()
454 {
455   // makes data from Digits
456
457   GetDigitsData(0)->Fill(fDigitsArray->GetEntriesFast()) ; 
458   TIter next(fDigitsArray) ; 
459   AliVZEROdigit *VZERODigit ; 
460   while ( (VZERODigit = dynamic_cast<AliVZEROdigit *>(next())) ) {
461     Int_t   PMNumber  = VZERODigit->PMNumber();         
462     GetDigitsData(PMNumber +1)->Fill( VZERODigit->Time()) ;    // in 100 of picoseconds
463     GetDigitsData(PMNumber +1+64)->Fill( VZERODigit->ADC()) ;
464   }  
465 }
466
467
468 //____________________________________________________________________________
469 void AliVZEROQADataMakerRec::MakeDigits(TTree *digitTree)
470 {
471   // makes data from Digit Tree
472         
473   if ( fDigitsArray ) 
474     fDigitsArray->Clear() ; 
475   else 
476     fDigitsArray = new TClonesArray("AliVZEROdigit", 1000) ; 
477   
478   TBranch * branch = digitTree->GetBranch("VZERODigit") ;
479   if ( ! branch ) {
480     AliWarning("VZERO branch in Digit Tree not found") ; 
481   } else {
482     branch->SetAddress(&fDigitsArray) ;
483     branch->GetEntry(0) ; 
484     MakeDigits() ; 
485   }  
486 }
487
488
489 //____________________________________________________________________________
490 void AliVZEROQADataMakerRec::MakeESDs(AliESDEvent * esd)
491 {
492   // Creates QA data from ESDs
493   
494   UInt_t eventType = esd->GetEventType();
495
496   switch (eventType){
497         case PHYSICS_EVENT:
498         AliESDVZERO *esdVZERO=esd->GetVZEROData();
499    
500         if (!esdVZERO) break;
501                   
502         GetESDsData(kCellMultiV0A)->Fill(esdVZERO->GetNbPMV0A());
503         GetESDsData(kCellMultiV0C)->Fill(esdVZERO->GetNbPMV0C());  
504         GetESDsData(kMIPMultiV0A)->Fill(esdVZERO->GetMTotV0A());
505         GetESDsData(kMIPMultiV0C)->Fill(esdVZERO->GetMTotV0C());  
506         
507         Float_t  timeV0A = 0., timeV0C = 0., diffTime;
508         Int_t   iTimeV0A = 0, iTimeV0C = 0;
509                 
510         for(Int_t i=0;i<64;i++) {
511                         GetESDsData(kMIPMultiChannel)->Fill((Float_t) i,(Float_t) esdVZERO->GetMultiplicity(i));
512                         GetESDsData(kChargeChannel)->Fill((Float_t) i,(Float_t) esdVZERO->GetAdc(i));
513                         if(esdVZERO->GetBBFlag(i)) GetESDsData(kBBFlag)->Fill((Float_t) i);
514                         if(esdVZERO->GetBGFlag(i)) GetESDsData(kBGFlag)->Fill((Float_t) i);
515
516                         Float_t time = (Float_t) esdVZERO->GetTime(i)/10.; //Convert in ns:  1 TDC channel = 100ps 
517                         GetESDsData(kTimeChannel)->Fill((Float_t) i,time);
518
519                         if(time>0.){
520                                 if (i<32) {
521                                         iTimeV0C++;
522                                         timeV0C += time;
523                                 }else{
524                                         iTimeV0A++;
525                                         timeV0A += time;
526                                 }
527                         }
528         }
529         if(iTimeV0A>0) timeV0A /= iTimeV0A; 
530         else timeV0A = -1.;
531         if(iTimeV0C>0) timeV0C /= iTimeV0C;
532         else timeV0C = -1.;
533         if(timeV0A<0. || timeV0C<0.) diffTime = -10000.;
534         else diffTime = timeV0A - timeV0C;
535                                 
536         GetESDsData(kESDV0ATime)->Fill(timeV0A);
537         GetESDsData(kESDV0CTime)->Fill(timeV0C);
538         GetESDsData(kESDDiffTime)->Fill(diffTime);
539                 
540         break;
541         }  
542   
543 }
544
545 //____________________________________________________________________________
546  void AliVZEROQADataMakerRec::MakeRaws(AliRawReader* rawReader)
547  {
548   // Fills histograms with Raws, computes average ADC values dynamically (pedestal subtracted)
549                   
550    // Check id histograms already created for this Event Specie
551    if ( ! GetRawsData(kPedestalInt0) )
552      InitRaws() ;
553
554    rawReader->Reset() ; 
555   AliVZERORawStream* rawStream  = new AliVZERORawStream(rawReader); 
556   rawStream->Next();
557   
558   eventTypeType eventType = rawReader->GetType();
559
560   Int_t    mulV0A = 0 ; 
561   Int_t    mulV0C = 0 ; 
562   Double_t timeV0A =0., timeV0C = 0.;
563   UInt_t   itimeV0A=0, itimeV0C=0;
564   Double_t chargeV0A=0., chargeV0C=0.;
565   Double_t mipV0A=0., mipV0C=0.;
566
567   Double_t diffTime=-100000.;
568
569   
570   switch (eventType){
571        case PHYSICS_EVENT:
572        Int_t  iFlag=0;
573        Int_t  pedestal;
574        Int_t  integrator;
575        Bool_t BBFlag;    
576        Bool_t BGFlag;    
577        UInt_t time, width;
578        Int_t  MBCharge, charge;
579        Int_t  offlineCh;
580        TH1D * hProj;
581
582        for(Int_t iChannel=0; iChannel<64; iChannel++) { // BEGIN : Loop over channels
583                    
584            offlineCh = rawStream->GetOfflineChannel(iChannel);
585                    
586            // Fill Pedestal histograms
587            
588            for(Int_t j=15; j<21; j++) {
589                        if((rawStream->GetBGFlag(iChannel,j) || rawStream->GetBBFlag(iChannel,j))) iFlag++;
590            }
591
592            if(iFlag == 0){ //No Flag found
593                        for(Int_t j=15; j<21; j++){
594                                pedestal=rawStream->GetPedestal(iChannel, j);
595                                integrator = rawStream->GetIntegratorFlag(iChannel, j);
596
597                                GetRawsData((integrator == 0 ? kPedestalInt0 : kPedestalInt1))->Fill(offlineCh,pedestal);
598                                GetRawsData((integrator == 0 ? kPedestalCycleInt0 : kPedestalCycleInt1))->Fill(offlineCh,pedestal);
599                        }
600             }
601
602            // Fill Charge EoI histograms
603            
604            // Look for the maximum in the LHC clock train
605            charge = 0;
606            Int_t iClock  = 0;
607            Int_t iCharge = 0;
608            for(Int_t iEvent=0; iEvent<21; iEvent++){
609                iCharge = rawStream->GetPedestal(iChannel,iEvent);
610                if(iCharge>charge)  {
611                        charge = iCharge;
612                        iClock = iEvent;
613                    }
614            }   // End of maximum searching procedure
615
616            integrator    = rawStream->GetIntegratorFlag(iChannel,iClock);
617            BBFlag        = rawStream->GetBBFlag(iChannel, iClock);
618            BGFlag        = rawStream->GetBGFlag(iChannel,iClock );
619
620            GetRawsData((integrator == 0 ? kChargeEoIInt0 : kChargeEoIInt1))->Fill(offlineCh,charge);
621            if(BBFlag) GetRawsData((integrator == 0 ? kChargeEoIBBInt0 : kChargeEoIBBInt1))->Fill(offlineCh,charge);
622            if(BGFlag) GetRawsData((integrator == 0 ? kChargeEoIBGInt0 : kChargeEoIBGInt1))->Fill(offlineCh,charge);
623
624            hProj = ((TH2I*)GetRawsData((integrator == 0 ? kPedestalInt0 : kPedestalInt1)))->ProjectionY("",offlineCh+1,offlineCh+1);
625            Double_t ped   = hProj->GetMean();
626            Double_t sigma = hProj->GetRMS();
627            delete hProj;
628
629            Double_t chargeEoI = charge - ped;
630                    
631            // Calculation of the number of MIP
632            Double_t mipEoI = chargeEoI * fCalibData->GetMIPperADC(offlineCh);
633
634                    
635            if(charge<1023 && chargeEoI > 5.*sigma){ 
636                    ((TH2I*)GetRawsData((integrator == 0 ? kChargeEoICycleInt0 : kChargeEoICycleInt1)))->Fill(offlineCh,chargeEoI);
637                    ((TH2D*)GetRawsData(kRawMIPChannel))->Fill(offlineCh,mipEoI);
638                    if(offlineCh<32) {
639                                    mulV0C++;
640                                    chargeV0C += chargeEoI;
641                                    mipV0C += mipEoI;
642                    } else {
643                                    mulV0A++;
644                                    chargeV0A += chargeEoI;
645                                    mipV0A += mipEoI;
646                    }
647            }
648
649            // Fill Charge Minimum Bias Histograms
650                    
651            int idx;
652            for(Int_t iBunch=0; iBunch<10; iBunch++){
653                            integrator = rawStream->GetIntMBFlag(iChannel, iBunch);
654                            BBFlag     = rawStream->GetBBMBFlag(iChannel, iBunch);
655                            BGFlag     = rawStream->GetBGMBFlag(iChannel, iBunch);
656                            MBCharge   = rawStream->GetChargeMB(iChannel, iBunch);
657
658                            if(integrator==0){
659                                    if(BBFlag==0){
660                                            if(BGFlag==0) idx = kChargeMBBB0BG0Int0;
661                                            else idx = kChargeMBBB0BG1Int0;
662                                    } else {
663                                            if(BGFlag==0) idx = kChargeMBBB1BG0Int0;
664                                            else idx = kChargeMBBB1BG1Int0;
665                                    }
666                            } else {
667                                    if(BBFlag==0){
668                                            if(BGFlag==0) idx = kChargeMBBB0BG0Int1;
669                                            else idx = kChargeMBBB0BG1Int1;
670                                    } else {
671                                            if(BGFlag==0) idx = kChargeMBBB1BG0Int1;
672                                            else idx = kChargeMBBB1BG1Int1;
673                                    }
674                            }
675                            GetRawsData(idx)->Fill(offlineCh,MBCharge);
676           }   
677
678           // Fill HPTDC Time Histograms
679
680            BBFlag   = rawStream->GetBBFlag(iChannel, 10);
681            BGFlag   = rawStream->GetBGFlag(iChannel, 10);
682            time     = rawStream->GetTime(iChannel);
683            width    = rawStream->GetWidth(iChannel);
684
685            if(time>0.){
686                       if (offlineCh<32) {
687                                    itimeV0C++;
688                                    timeV0C += time;
689                       }else{
690                                    itimeV0A++;
691                                    timeV0A += time;
692                       }
693            }
694            GetRawsData(kHPTDCTime)->Fill(offlineCh,time);
695            GetRawsData(kWidth)->Fill(offlineCh,width);
696            if(BBFlag) {
697                   GetRawsData(kHPTDCTimeBB)->Fill(offlineCh,time);
698                   GetRawsData(kWidthBB)->Fill(offlineCh,width);
699            }
700            if(BGFlag) {
701                   GetRawsData(kHPTDCTimeBG)->Fill(offlineCh,time);
702                   GetRawsData(kWidthBG)->Fill(offlineCh,width);
703            }
704
705            // Fill Flag and Charge Versus LHC-Clock histograms
706            
707            for(Int_t iEvent=0; iEvent<21; iEvent++){
708                charge = rawStream->GetPedestal(iChannel,iEvent);
709                integrator = rawStream->GetIntegratorFlag(iChannel,iEvent);
710                BBFlag     = rawStream->GetBBFlag(iChannel, iEvent);
711                BGFlag     = rawStream->GetBGFlag(iChannel,iEvent );
712
713                ((TH2*) GetRawsData((integrator == 0 ? kChargeVsClockInt0 : kChargeVsClockInt1 )))->Fill(offlineCh,(float)iEvent-10,(float)charge);
714                ((TH2*) GetRawsData(kBBFlagVsClock))->Fill(offlineCh,(float)iEvent-10,(float)BBFlag);
715                ((TH2*) GetRawsData(kBGFlagVsClock))->Fill(offlineCh,(float)iEvent-10,(float)BGFlag);
716            }
717
718        }// END of Loop over channels
719
720             if(itimeV0A>0) timeV0A /= (itimeV0A * 10); // itimeV0A Channels and divide by 10 to have the result in ns because 1 TDC Channel = 100 ps
721             else timeV0A = -1.;
722             if(itimeV0C>0) timeV0C /= (itimeV0C * 10);
723             else timeV0C = -1.;
724             if(timeV0A<0. || timeV0C<0.) diffTime = -10000.;
725             else diffTime = timeV0A - timeV0C;
726                 
727             GetRawsData(kV0ATime)->Fill(timeV0A);
728             GetRawsData(kV0CTime)->Fill(timeV0C);
729             GetRawsData(kDiffTime)->Fill(diffTime);
730                 
731             GetRawsData(kMultiV0A)->Fill(mulV0A);
732             GetRawsData(kMultiV0C)->Fill(mulV0C);
733
734             GetRawsData(kChargeV0A)->Fill(chargeV0A);
735             GetRawsData(kChargeV0C)->Fill(chargeV0C);
736             GetRawsData(kChargeV0)->Fill(chargeV0A + chargeV0C);
737                 
738             GetRawsData(kRawMIPV0A)->Fill(mipV0A);
739             GetRawsData(kRawMIPV0C)->Fill(mipV0C);
740             GetRawsData(kRawMIPV0)->Fill(mipV0A + mipV0C);
741             break;
742             
743         } // END of SWITCH : EVENT TYPE 
744         
745         fEvent++; 
746         TParameter<double> * p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kMultiV0A)->GetName()))) ; 
747         if (p) p->SetVal((double)mulV0A) ; 
748
749         p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kMultiV0C)->GetName()))) ; 
750         if (p) p->SetVal((double)mulV0C) ;                     
751
752         p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kChargeV0A)->GetName()))) ; 
753         if (p) p->SetVal((double)chargeV0A) ; 
754
755         p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kChargeV0C)->GetName()))) ; 
756         if (p) p->SetVal((double)chargeV0C) ;                     
757
758         p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kChargeV0)->GetName()))) ; 
759         if (p) p->SetVal((double)(chargeV0A + chargeV0C)) ;                     
760         
761         p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kRawMIPV0A)->GetName()))) ; 
762         if (p) p->SetVal((double)mipV0A) ; 
763         
764         p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kRawMIPV0C)->GetName()))) ; 
765         if (p) p->SetVal((double)mipV0C) ;                     
766         
767         p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kRawMIPV0)->GetName()))) ; 
768         if (p) p->SetVal((double)(mipV0A + mipV0C)) ;                     
769         
770         p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kV0ATime)->GetName()))) ; 
771         if (p) p->SetVal((double)timeV0A) ; 
772         
773         p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kV0CTime)->GetName()))) ; 
774         if (p) p->SetVal((double)timeV0C) ;                     
775         
776         p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(AliQAv1::kRAWS).Data(), GetRawsData(kDiffTime)->GetName()))) ; 
777         if (p) p->SetVal((double)diffTime) ;                     
778         
779         delete rawStream; rawStream = 0x0;      
780
781
782  }
783
784 //____________________________________________________________________________ 
785 void AliVZEROQADataMakerRec::StartOfDetectorCycle()
786 {
787   // Detector specific actions at start of cycle
788   
789   // Reset of the histogram used - to have the trend versus time -
790  
791   fCalibData = GetCalibData();
792   
793   TH1* h;
794   h = GetRawsData(kPedestalCycleInt0);
795   if(h) h->Reset();
796   h = GetRawsData(kPedestalCycleInt1); 
797   if(h) h->Reset();
798   h = GetRawsData(kChargeEoICycleInt0);
799   if(h) h->Reset();
800   h = GetRawsData(kChargeEoICycleInt1);
801   if(h) h->Reset();
802
803 }