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