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