]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCdataQA.cxx
Add three particle correlation loop for correction of auto-correlaton bias
[u/mrichter/AliRoot.git] / TPC / AliTPCdataQA.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 /* $Id$ */
18
19 /*
20   July 2011:
21
22   Changes to accomodate updates of general DQM/QA changes to have per trigger
23   histograms (for a given event specie).
24
25   AliTPCdataQA has a new flag for only keeping DQM info event by
26   event!
27   The expert/DA functionality has been kept exactly the same!
28   
29
30   June 2010
31
32   This update should solve two problems mainly:
33   * The vs event histograms have been limited to a fixed size for the
34   DQM. The 500k seemed to be a big size but is no longer so, so we
35   need to dynamically expand the range. The non-trivial point is that
36   we also have to do it for the copy owned by AliTPCQADataMakerRec.
37   * The amoreGui now remembers the x-range of the first visualization so
38   the trick of setting the relevant event range as the histogram is
39   filled no longer works.
40
41   The fix is a bit crude but avoids creating a new histogram. Instead
42   the range is expanded (max events and events per bin is doubled) but
43   the number of bins is kept constant! In this way we can change just
44   double the max of the X-axis of the hist and rebin the data. The
45   same can easily be done for the copy owned by AliTPCQADataMakerRec.
46
47   CAUTION:
48   If we change the number of bins we could crash the whole system
49   because ROOT does not create space for extra bins! (but we do not do
50   this). In that way it is a crude solution.
51   The rebinning in the code only works for an even number of bins.
52
53   In addition to the above a bug in the reading of the config file was
54   also found and corrected. fAdcMax was set instead of fEventsPerBin.
55
56   Finally cout was changes to AliInfo.
57
58   February 2008
59
60   The code has been heavily modified so that now the RAW data is
61   "expanded" for each sector and stored in a big signal array. Then a
62   simple version of the code in AliTPCclustererMI is used to identify
63   the local maxima and these are then used for QA. This gives a better
64   estimate of the charge (both max and total) and also limits the
65   effect of noise.
66
67   Implementation:
68
69   In Update the RAW signals >= 3 ADC channels are stored in the arrays.
70   
71   There are 3 arrays:
72   Float_t** fAllBins       2d array [row][bin(pad, time)] ADC signal
73   Int_t**   fAllSigBins    2d array [row][signal#] bin(with signal)
74   Int_t*    fAllNSigBins;  1d array [row] Nsignals
75
76   This is done sector by sector.
77
78   When data from a new sector is encountered, the method
79   FindLocalMaxima is called on the data from the previous sector, and
80   the calibration/data objects are updated with the "cluster"
81   info. Finally the arrays are cleared.
82
83   The requirements for a local maxima is:
84   Charge in bin is >= 5 ADC channels.
85   Charge in bin is larger than all the 8 neighboring bins.
86   At least one of the two pad neighbors has a signal.
87   At least one of the two time neighbors has a signal.
88
89   Before accessing the data it is expected that the Analyse method is
90   called. This normalizes some of the data objects to per event or per
91   cluster. 
92   If more data is passed to the class after Analyse has been called
93   the normalization is reversed and Analyse has to be called again.
94 */
95
96
97 //Root includes
98 #include <TH1F.h>
99 #include <TString.h>
100 #include <TMath.h>
101 #include <TDirectory.h>
102 #include <TFile.h>
103 #include <TError.h>
104 #include <TMap.h>
105 #include <TProfile.h>
106 //AliRoot includes
107 #include "AliRawReader.h"
108 #include "AliRawReaderRoot.h"
109 #include "AliRawReaderDate.h"
110 #include "AliTPCRawStream.h"
111 #include "AliTPCRawStreamV3.h"
112 #include "AliTPCCalROC.h"
113 #include "AliTPCROC.h"
114 #include "AliMathBase.h"
115 #include "TTreeStream.h"
116
117 //date
118 #include "event.h"
119 #include "AliTPCCalPad.h"
120 #include "AliTPCPreprocessorOnline.h"
121
122 //header file
123 #include "AliTPCdataQA.h"
124 #include "AliLog.h"
125
126
127 ClassImp(AliTPCdataQA)
128
129 AliTPCdataQA::AliTPCdataQA() : /*FOLD00*/  
130   fFirstTimeBin(60),
131   fLastTimeBin(1000),
132   fAdcMin(1),
133   fAdcMax(100),
134   fMapping(NULL),
135   fPedestal(0),
136   fNoise(0),
137   fNLocalMaxima(0),
138   fMaxCharge(0),
139   fMeanCharge(0),
140   fNoThreshold(0),
141   fNTimeBins(0),
142   fNPads(0),
143   fTimePosition(0),
144   fOverThreshold10(0),
145   fOverThreshold20(0),
146   fOverThreshold30(0),
147   fHistQVsTimeSideA(0),
148   fHistQVsTimeSideC(0),
149   fHistQMaxVsTimeSideA(0),
150   fHistQMaxVsTimeSideC(0),
151   fHistOccupancyVsEvent(0),
152   fHistNclustersVsEvent(0),
153   fEventCounter(0),
154   fIsAnalysed(kFALSE),
155   fMaxEvents(500000),           // Max events for event histograms
156   fEventsPerBin(1000),          // Events per bin for event histograms
157   fSignalCounter(0),            // Signal counter
158   fClusterCounter(0),           // Cluster counter
159   fAllBins(0),
160   fAllSigBins(0),
161   fAllNSigBins(0),
162   fRowsMax(0),
163   fPadsMax(0),
164   fTimeBinsMax(0),
165   fIsDQM(kFALSE),
166   fHistOccVsSector(0x0),
167   fHistOcc2dVsSector(0x0),
168   fHistQVsSector(0x0),
169   fHistQmaxVsSector(0x0),
170   fOccVec(0x0),
171   fOccMaxVec(0x0),
172   fOccVecFine(0x0),
173   fOccMaxVecFine(0x0)
174 {
175   //
176   // default constructor
177   //
178 }
179
180 //_____________________________________________________________________
181 AliTPCdataQA::AliTPCdataQA(const AliTPCdataQA &ped) : /*FOLD00*/
182   TH1F(ped),
183   fFirstTimeBin(ped.GetFirstTimeBin()),
184   fLastTimeBin(ped.GetLastTimeBin()),
185   fAdcMin(ped.GetAdcMin()),
186   fAdcMax(ped.GetAdcMax()),
187   fMapping(NULL),
188   fPedestal(0),
189   fNoise(0),
190   fNLocalMaxima(0),
191   fMaxCharge(0),
192   fMeanCharge(0),
193   fNoThreshold(0),
194   fNTimeBins(0),
195   fNPads(0),
196   fTimePosition(0),
197   fOverThreshold10(0),
198   fOverThreshold20(0),
199   fOverThreshold30(0),
200   fHistQVsTimeSideA(0),
201   fHistQVsTimeSideC(0),
202   fHistQMaxVsTimeSideA(0),
203   fHistQMaxVsTimeSideC(0),
204   fHistOccupancyVsEvent(0),
205   fHistNclustersVsEvent(0),
206   fEventCounter(ped.GetEventCounter()),
207   fIsAnalysed(ped.GetIsAnalysed()),
208   fMaxEvents(ped.GetMaxEvents()),
209   fEventsPerBin(ped.GetEventsPerBin()),
210   fSignalCounter(ped.GetSignalCounter()),
211   fClusterCounter(ped.GetClusterCounter()),
212   fAllBins(0),
213   fAllSigBins(0),
214   fAllNSigBins(0),
215   fRowsMax(0),
216   fPadsMax(0),
217   fTimeBinsMax(0),
218   fIsDQM(ped.GetIsDQM()),
219   fHistOccVsSector(0x0),
220   fHistOcc2dVsSector(0x0),
221   fHistQVsSector(0x0),
222   fHistQmaxVsSector(0x0),
223   fOccVec(0x0),
224   fOccMaxVec(0x0),
225   fOccVecFine(0x0),
226   fOccMaxVecFine(0x0)
227 {
228   //
229   // copy constructor
230   //
231   if(ped.GetNLocalMaxima())
232     fNLocalMaxima  = new AliTPCCalPad(*ped.GetNLocalMaxima());
233   if(ped.GetMaxCharge())
234     fMaxCharge      = new AliTPCCalPad(*ped.GetMaxCharge());
235   if(ped.GetMeanCharge())
236     fMeanCharge     = new AliTPCCalPad(*ped.GetMeanCharge());
237   if(ped.GetNoThreshold())
238     fNoThreshold  = new AliTPCCalPad(*ped.GetNoThreshold());
239   if(ped.GetNTimeBins())
240     fNTimeBins  = new AliTPCCalPad(*ped.GetNTimeBins());
241   if(ped.GetNPads())
242     fNPads  = new AliTPCCalPad(*ped.GetNPads());
243   if(ped.GetTimePosition())
244     fTimePosition  = new AliTPCCalPad(*ped.GetTimePosition());
245   if(ped.GetOverThreshold10())
246     fOverThreshold10  = new AliTPCCalPad(*ped.GetOverThreshold10());
247   if(ped.GetOverThreshold20())
248     fOverThreshold20  = new AliTPCCalPad(*ped.GetOverThreshold20());
249   if(ped.GetOverThreshold30())
250     fOverThreshold30  = new AliTPCCalPad(*ped.GetOverThreshold30());
251   if(ped.GetHistQVsTimeSideA()) {
252     fHistQVsTimeSideA = new TProfile(*ped.GetHistQVsTimeSideA());
253     fHistQVsTimeSideA->SetDirectory(0);
254   }
255   if(ped.GetHistQVsTimeSideC()) {
256     fHistQVsTimeSideC = new TProfile(*ped.GetHistQVsTimeSideC());
257     fHistQVsTimeSideC->SetDirectory(0);
258   }
259   if(ped.GetHistQMaxVsTimeSideA()) {
260     fHistQMaxVsTimeSideA = new TProfile(*ped.GetHistQMaxVsTimeSideA());
261     fHistQMaxVsTimeSideA->SetDirectory(0);
262   }
263   if(ped.GetHistQMaxVsTimeSideC()) {
264     fHistQMaxVsTimeSideC = new TProfile(*ped.GetHistQMaxVsTimeSideC());
265     fHistQMaxVsTimeSideC->SetDirectory(0);
266   }
267   if(ped.GetHistOccupancyVsEventConst()) {
268     fHistOccupancyVsEvent  = new TH1F(*ped.GetHistOccupancyVsEventConst());
269     fHistOccupancyVsEvent->SetDirectory(0);
270   }
271   if(ped.GetHistNclustersVsEventConst()) {
272     fHistNclustersVsEvent  = new TH1F(*ped.GetHistNclustersVsEventConst());
273     fHistNclustersVsEvent->SetDirectory(0);
274   }
275 }
276
277 //_____________________________________________________________________
278 AliTPCdataQA::AliTPCdataQA(const TMap *config) : /*FOLD00*/  
279   TH1F("TPCRAW","TPCRAW",100,0,100),
280   fFirstTimeBin(60),
281   fLastTimeBin(1000),
282   fAdcMin(1),
283   fAdcMax(100),
284   fMapping(NULL),
285   fPedestal(0),
286   fNoise(0),
287   fNLocalMaxima(0),
288   fMaxCharge(0),
289   fMeanCharge(0),
290   fNoThreshold(0),
291   fNTimeBins(0),
292   fNPads(0),
293   fTimePosition(0),
294   fOverThreshold10(0),
295   fOverThreshold20(0),
296   fOverThreshold30(0),
297   fHistQVsTimeSideA(0),
298   fHistQVsTimeSideC(0),
299   fHistQMaxVsTimeSideA(0),
300   fHistQMaxVsTimeSideC(0),
301   fHistOccupancyVsEvent(0),
302   fHistNclustersVsEvent(0),
303   fEventCounter(0),
304   fIsAnalysed(kFALSE),
305   fMaxEvents(500000),
306   fEventsPerBin(1000),
307   fSignalCounter(0),
308   fClusterCounter(0),
309   fAllBins(0),
310   fAllSigBins(0),
311   fAllNSigBins(0),
312   fRowsMax(0),
313   fPadsMax(0),
314   fTimeBinsMax(0),
315   fIsDQM(kFALSE),
316   fHistOccVsSector(0x0),
317   fHistOcc2dVsSector(0x0),
318   fHistQVsSector(0x0),
319   fHistQmaxVsSector(0x0),
320   fOccVec(0x0),
321   fOccMaxVec(0x0),
322   fOccVecFine(0x0),
323   fOccMaxVecFine(0x0)
324 {
325   //
326   // default constructor
327   //
328   if (config->GetValue("FirstTimeBin")) fFirstTimeBin = ((TObjString*)config->GetValue("FirstTimeBin"))->GetString().Atoi();
329   if (config->GetValue("LastTimeBin"))  fLastTimeBin = ((TObjString*)config->GetValue("LastTimeBin"))->GetString().Atoi();
330   if (config->GetValue("AdcMin"))       fAdcMin = ((TObjString*)config->GetValue("AdcMin"))->GetString().Atoi();
331   if (config->GetValue("AdcMax"))       fAdcMax = ((TObjString*)config->GetValue("AdcMax"))->GetString().Atoi();
332   if (config->GetValue("MaxEvents"))    fMaxEvents = ((TObjString*)config->GetValue("MaxEvents"))->GetString().Atoi();
333   if (config->GetValue("EventsPerBin")) fEventsPerBin = ((TObjString*)config->GetValue("EventsPerBin"))->GetString().Atoi();
334 }
335
336 //_____________________________________________________________________
337 AliTPCdataQA& AliTPCdataQA::operator = (const  AliTPCdataQA &source)
338 {
339   //
340   // assignment operator
341   //
342   if (&source == this) return *this;
343   new (this) AliTPCdataQA(source);
344
345   return *this;
346 }
347
348
349 //_____________________________________________________________________
350 AliTPCdataQA::~AliTPCdataQA() /*FOLD00*/
351 {
352   //
353   // destructor
354   //
355
356   // do not delete fMapping, because we do not own it.
357   // do not delete fMapping, because we do not own it.
358   // do not delete fNoise and fPedestal, because we do not own them.
359
360   delete fNLocalMaxima;
361   delete fMaxCharge;
362   delete fMeanCharge;
363   delete fNoThreshold;
364   delete fNTimeBins;
365   delete fNPads;
366   delete fTimePosition;
367   delete fOverThreshold10;
368   delete fOverThreshold20;
369   delete fOverThreshold30;
370   delete fHistQVsTimeSideA;
371   delete fHistQVsTimeSideC;
372   delete fHistQMaxVsTimeSideA;
373   delete fHistQMaxVsTimeSideC;
374   delete fHistOccupancyVsEvent;
375   delete fHistNclustersVsEvent;
376
377   // DQM
378   delete fHistOccVsSector;
379   delete fHistOcc2dVsSector;
380   delete fHistQVsSector;
381   delete fHistQmaxVsSector;
382   delete fOccVec;
383   delete fOccMaxVec;
384   delete fOccVecFine;
385   delete fOccMaxVecFine;
386   
387   for (Int_t iRow = 0; iRow < fRowsMax; iRow++) {
388     delete [] fAllBins[iRow];
389     delete [] fAllSigBins[iRow];
390   }  
391   delete [] fAllBins;
392   delete [] fAllSigBins;
393   delete [] fAllNSigBins;
394 }
395
396 //_____________________________________________________________________
397 TH1F* AliTPCdataQA::GetHistOccupancyVsEvent()
398 {
399   //
400   // Create Occupancy vs event histogram
401   // (we create this histogram differently then the other histograms
402   //  because this we want to be able to access and copy
403   //  from AliTPCQAMakerRec before it normally would be created)
404   //
405   if(!fHistOccupancyVsEvent) {
406
407     Int_t nBins = fMaxEvents/fEventsPerBin;
408     fHistOccupancyVsEvent = new TH1F("hOccupancyVsEvent", "Occupancy vs event number (~time); Event number; Occupancy", nBins, 0, nBins*fEventsPerBin);
409     fHistOccupancyVsEvent->SetDirectory(0);
410   }
411   
412   return fHistOccupancyVsEvent;
413 }
414
415 //_____________________________________________________________________
416 TH1F* AliTPCdataQA::GetHistNclustersVsEvent()
417 {
418   //
419   // Create Nclusters vs event histogram
420   // (we create this histogram differently then the other histograms
421   //  because this we want to be able to access and copy
422   //  from AliTPCQAMakerRec before it normally would be created)
423   //
424   if(!fHistNclustersVsEvent) {
425
426     Int_t nBins = fMaxEvents/fEventsPerBin;
427     fHistNclustersVsEvent = new TH1F("hNclustersVsEvent", "Nclusters vs event number (~time); Event number; Nclusters per event", nBins, 0, nBins*fEventsPerBin);
428     fHistNclustersVsEvent->SetDirectory(0);
429   }
430   
431   return fHistNclustersVsEvent;
432 }
433
434 //_____________________________________________________________________
435 void AliTPCdataQA::UpdateEventHistograms()
436 {
437   // Update histograms that display occupancy and 
438   // number of clusters as a function of number of 
439   // events
440   if (!fHistOccupancyVsEvent)
441     GetHistOccupancyVsEvent();
442   if (!fHistNclustersVsEvent)
443     GetHistNclustersVsEvent();
444   
445   if(fEventCounter > fMaxEvents) {
446     
447     // we have to expand the histogram to handle the larger number of
448     // events. The way it is done now is to double the range and the
449     // number of events per bin (so the number of histogram bins stays
450     // constant)
451     fEventsPerBin *= 2;
452     fMaxEvents *= 2;
453
454     // Change histogram limits
455     const Int_t nBins = fHistOccupancyVsEvent->GetXaxis()->GetNbins();
456     fHistOccupancyVsEvent->GetXaxis()->Set(nBins, fHistOccupancyVsEvent->GetXaxis()->GetNbins(), fMaxEvents);
457     fHistNclustersVsEvent->GetXaxis()->Set(nBins, fHistNclustersVsEvent->GetXaxis()->GetNbins(), fMaxEvents);
458
459     // Rebin the histogram
460     for(Int_t bin = 1; bin <= nBins; bin+=2) {
461
462       Int_t newBin = TMath::Nint(Float_t(bin+1)/2.0);
463       Float_t newContent = (fHistOccupancyVsEvent->GetBinContent(bin)+
464                             fHistOccupancyVsEvent->GetBinContent(bin+1))/2.0;
465       fHistOccupancyVsEvent->SetBinContent(newBin, newContent); 
466
467       newContent = (fHistNclustersVsEvent->GetBinContent(bin)+
468                     fHistNclustersVsEvent->GetBinContent(bin+1))/2.0;
469       fHistNclustersVsEvent->SetBinContent(newBin, newContent); 
470     }
471
472     // Set the remaining bins to 0
473     Int_t lastHalf = nBins/2 +1;
474     for(Int_t bin = lastHalf; bin <= nBins; bin++) {
475
476       fHistOccupancyVsEvent->SetBinContent(bin, 0); 
477       fHistNclustersVsEvent->SetBinContent(bin, 0); 
478     }
479
480     // In this case we should nut update but wait untill the new
481     // number of events per bin is reached!
482     return;
483   }
484
485   const Int_t bin = TMath::Nint(Float_t(fEventCounter)/fEventsPerBin);
486
487   Float_t averageOccupancy =
488     Float_t(fSignalCounter)/fEventsPerBin/(fLastTimeBin - fFirstTimeBin +1.0)
489     / 570132.0; // 570,132 is number of pads
490   fHistOccupancyVsEvent->SetBinContent(bin, averageOccupancy);
491   fSignalCounter = 0;
492   
493   Float_t averageNclusters =
494     Float_t(fClusterCounter)/fEventsPerBin;
495   fHistNclustersVsEvent->SetBinContent(bin, averageNclusters);
496   fClusterCounter = 0;
497 }
498
499 //_____________________________________________________________________
500 Bool_t AliTPCdataQA::ProcessEvent(AliTPCRawStreamV3 *const rawStreamV3)
501 {
502   //
503   // Event Processing loop - AliTPCRawStreamV3
504   //
505   Bool_t withInput = kFALSE;
506   Int_t nSignals = 0;
507   Int_t lastSector = -1;
508   
509   while ( rawStreamV3->NextDDL() ){
510
511     while ( rawStreamV3->NextChannel() ){
512     
513       Int_t iSector = rawStreamV3->GetSector(); //  current sector
514       Int_t iRow    = rawStreamV3->GetRow();    //  current row
515       Int_t iPad    = rawStreamV3->GetPad();    //  current pad
516       Int_t iPatch  = rawStreamV3->GetPatchIndex(); //  current patch
517       Int_t iBranch = rawStreamV3->GetBranch();    //  current branch
518       if (iRow<0 || iPad<0) continue;
519       // Call local maxima finder if the data is in a new sector
520       if(iSector != lastSector) {
521         
522         if(nSignals>0)
523           FindLocalMaxima(lastSector);
524         
525         CleanArrays();
526         lastSector = iSector;
527         nSignals = 0;
528       }
529       
530       while ( rawStreamV3->NextBunch() ){
531
532         Int_t  startTbin    = (Int_t)rawStreamV3->GetStartTimeBin();
533         Int_t  bunchlength  = (Int_t)rawStreamV3->GetBunchLength();
534         const UShort_t *sig = rawStreamV3->GetSignals();
535         
536         for (Int_t iTimeBin = 0; iTimeBin<bunchlength; iTimeBin++){
537           Float_t signal=(Float_t)sig[iTimeBin];
538           nSignals += Update(iSector,iRow,iPad,startTbin--,signal, iPatch, iBranch);
539           withInput = kTRUE;
540         }
541       }
542     }
543   }
544     
545   if (lastSector>=0&&nSignals>0)  
546     FindLocalMaxima(lastSector);
547   
548   return withInput;
549 }
550
551 //_____________________________________________________________________
552 Bool_t AliTPCdataQA::ProcessEvent(AliRawReader *const rawReader)
553 {
554   //
555   //  Event processing loop - AliRawReader
556   //
557   AliTPCRawStreamV3 rawStreamV3(rawReader,(AliAltroMapping**)fMapping);
558   Bool_t res=ProcessEvent(&rawStreamV3);
559   if(res) {
560     fEventCounter++; // only increment event counter if there is TPC data
561
562     if(fEventCounter%fEventsPerBin==0) 
563       UpdateEventHistograms();
564   }
565   return res;
566 }
567
568 //_____________________________________________________________________
569 Bool_t AliTPCdataQA::ProcessEvent(AliTPCRawStream *const rawStream)
570 {
571   //
572   // Event Processing loop - AliTPCRawStream
573   //
574
575   Bool_t withInput = kFALSE;
576   Int_t nSignals = 0;
577   Int_t lastSector = -1;
578
579   while (rawStream->Next()) {
580
581     Int_t iSector  = rawStream->GetSector();      //  current ROC
582     Int_t iRow     = rawStream->GetRow();         //  current row
583     Int_t iPad     = rawStream->GetPad();         //  current pad
584     Int_t iTimeBin = rawStream->GetTime();        //  current time bin
585     Float_t signal = rawStream->GetSignal();      //  current ADC signal
586     
587     // Call local maxima finder if the data is in a new sector
588     if(iSector != lastSector) {
589       
590       if(nSignals>0)
591         FindLocalMaxima(lastSector);
592       
593       CleanArrays();
594       lastSector = iSector;
595       nSignals = 0;
596     }
597     
598     // Sometimes iRow==-1 if there is problems to read the data
599     if(iRow>=0) {
600       nSignals += Update(iSector,iRow,iPad,iTimeBin,signal);
601       withInput = kTRUE;
602     }
603   }
604
605   if (lastSector>=0&&nSignals>0)
606     FindLocalMaxima(lastSector);
607   
608   return withInput;
609 }
610
611
612 //_____________________________________________________________________
613 Bool_t AliTPCdataQA::ProcessEventOld(AliRawReader *const rawReader)
614 {
615   //
616   //  Event processing loop - AliRawReader
617   //
618
619   // if fMapping is NULL the rawstream will crate its own mapping
620   AliTPCRawStream rawStream(rawReader, (AliAltroMapping**)fMapping);
621   rawReader->Select("TPC");
622   Bool_t res =  ProcessEvent(&rawStream);
623
624   if(res)
625     fEventCounter++; // only increment event counter if there is TPC data
626                      // otherwise Analyse (called in QA) fails
627   return res;
628 }
629
630
631 //_____________________________________________________________________
632 Bool_t AliTPCdataQA::ProcessEvent(eventHeaderStruct *const event)
633 {
634   //
635   //  process date event
636   //
637
638   AliRawReaderDate rawReader((void*)event);
639   Bool_t result=ProcessEvent(&rawReader);
640   return result;
641 }
642
643
644
645 //_____________________________________________________________________
646 void AliTPCdataQA::DumpToFile(const Char_t *filename, const Char_t *dir, Bool_t append) /*FOLD00*/
647 {
648   //
649   //  Write class to file
650   //
651
652   TString sDir(dir);
653   TString option;
654
655   if ( append )
656     option = "update";
657   else
658     option = "recreate";
659
660   TDirectory *backup = gDirectory;
661   TFile f(filename,option.Data());
662   f.cd();
663   if ( !sDir.IsNull() ){
664     f.mkdir(sDir.Data());
665     f.cd(sDir);
666   }
667   this->Write();
668   f.Close();
669
670   if ( backup ) backup->cd();
671 }
672
673
674 //_____________________________________________________________________
675 Int_t AliTPCdataQA::Update(const Int_t iSector, /*FOLD00*/
676                            const Int_t iRow,
677                            const Int_t iPad,
678                            const Int_t iTimeBin,
679                            Float_t signal,
680                            const Int_t iPatch,
681                            const Int_t iBranch)
682 {
683   //
684   // Signal filling method
685   //
686   
687   //
688   // Define the calibration objects the first time Update is called
689   // NB! This has to be done first even if the data is rejected by the time 
690   // cut to make sure that the objects are available in Analyse
691   //
692   if(!fIsDQM) {
693
694     if (!fNLocalMaxima) fNLocalMaxima = new AliTPCCalPad("NLocalMaxima","NLocalMaxima");
695     if (!fMaxCharge) fMaxCharge = new AliTPCCalPad("MaxCharge","MaxCharge");
696     if (!fMeanCharge) fMeanCharge = new AliTPCCalPad("MeanCharge","MeanCharge");
697     if (!fNoThreshold) fNoThreshold = new AliTPCCalPad("NoThreshold","NoThreshold");
698     if (!fNTimeBins) fNTimeBins = new AliTPCCalPad("NTimeBins","NTimeBins");
699     if (!fNPads) fNPads = new AliTPCCalPad("NPads","NPads");
700     if (!fTimePosition) fTimePosition = new AliTPCCalPad("TimePosition","TimePosition");
701     if (!fOverThreshold10) fOverThreshold10 = new AliTPCCalPad("OverThreshold10","OverThreshold10");
702     if (!fOverThreshold20) fOverThreshold20 = new AliTPCCalPad("OverThreshold20","OverThreshold20");
703     if (!fOverThreshold30) fOverThreshold30 = new AliTPCCalPad("OverThreshold30","OverThreshold30");
704     if (!fHistQVsTimeSideA) {
705       fHistQVsTimeSideA  = new TProfile("hQVsTimeSideA", "Q vs time (side A); Time [Timebin]; Q [ADC ch]", 100, 0, 1000);
706       fHistQVsTimeSideA->SetDirectory(0);
707     }
708     if (!fHistQVsTimeSideC) {
709       fHistQVsTimeSideC  = new TProfile("hQVsTimeSideC", "Q vs time (side C); Time [Timebin]; Q [ADC ch]", 100, 0, 1000);
710       fHistQVsTimeSideC->SetDirectory(0);
711   }
712     if (!fHistQMaxVsTimeSideA) {
713       fHistQMaxVsTimeSideA  = new TProfile("hQMaxVsTimeSideA", "Q_{MAX} vs time (side A); Time [Timebin]; Q_{MAX} [ADC ch]", 100, 0, 1000);
714       fHistQMaxVsTimeSideA->SetDirectory(0);
715     }
716     if (!fHistQMaxVsTimeSideC) {
717       fHistQMaxVsTimeSideC  = new TProfile("hQMaxVsTimeSideC", "Q_{MAX} vs time (side C); Time [Timebin]; Q_{MAX} [ADC ch]", 100, 0, 1000);
718       fHistQMaxVsTimeSideC->SetDirectory(0);
719     }
720   } else { // DQM histograms and array
721     
722     if (!fHistOccVsSector) {
723       fHistOccVsSector  = new TProfile("hOccVsSector", "Occupancy vs sector; Sector; Occupancy", 72, 0, 72);
724       fHistOccVsSector->SetDirectory(0);
725
726       fHistOcc2dVsSector  = new TProfile2D("hOcc2dVsSector", "Occupancy vs sector and patch; Sector; Patch", 72, 0, 36, 6, 0, 6);
727       fHistOcc2dVsSector->SetDirectory(0);
728
729       fHistQVsSector  = new TProfile("hQVsSector", "Q vs sector; Sector; Q [ADC ch]", 72, 0, 72);
730       fHistQVsSector->SetDirectory(0);
731
732       fHistQmaxVsSector  = new TProfile("hQmaxVsSector", "Qmax vs sector; Sector; Qmax [ADC ch]", 72, 0, 72);
733       fHistQmaxVsSector->SetDirectory(0);
734
735       fOccVec = new TArrayD(72);
736       for(Int_t i = 0; i < 72; i++)
737         fOccVec->GetArray()[i] = 0;
738
739       fOccMaxVec = new TArrayD(72);
740       const Double_t nTimeBins = fLastTimeBin - fFirstTimeBin +1;
741       for(Int_t i = 0; i < 72; i++)
742         
743         if(i<36) // IROCs (5504 pads)
744           fOccMaxVec->GetArray()[i] = nTimeBins*5504;
745         else     // OROCs (9984 pads)
746           fOccMaxVec->GetArray()[i] = nTimeBins*9984;
747
748       // 12 branches for each full sector
749       const Int_t nBranches = 36*12;
750       fOccVecFine = new TArrayD(nBranches);
751       for(Int_t i = 0; i < nBranches; i++)
752         fOccVecFine->GetArray()[i] = 0;
753
754       // Pads per patch same for all sectors
755       Int_t nPads0[6] = {1152, 1536, 1152, 1280, 1280, 1280};
756       Int_t nPads1[6] = {1152, 1664, 1152, 1280, 1280, 1280};
757
758       fOccMaxVecFine = new TArrayD(nBranches);
759       for(Int_t i = 0; i < nBranches; i++) {
760         
761         const Int_t fullSector = Int_t(i/12);
762         Int_t branch = i - fullSector*12;
763         R__ASSERT(branch>=0 && branch<12);
764         
765         const Int_t patch = Int_t(branch/2);
766         branch -= patch*2;
767         
768         R__ASSERT(branch>=0 && branch<2);
769         if(branch == 0)
770           fOccMaxVecFine->GetArray()[i] = nTimeBins*nPads0[patch];
771         else     // OROCs (9984 pads)
772           fOccMaxVecFine->GetArray()[i] = nTimeBins*nPads1[patch];
773       }
774     }
775   }
776   // Make the arrays for expanding the data
777
778   if (!fAllBins)
779     MakeArrays();
780
781   //
782   // If Analyse has been previously called we need now to denormalize the data
783   // as more data is coming
784   //
785   if(fIsAnalysed == kTRUE && !fIsDQM) {
786     
787     const Int_t nTimeBins = fLastTimeBin - fFirstTimeBin +1;
788     const Float_t denormalization = Float_t(fEventCounter * nTimeBins);
789     fNoThreshold->Multiply(denormalization);  
790     
791     fMeanCharge->Multiply(fNLocalMaxima);
792     fMaxCharge->Multiply(fNLocalMaxima);
793     fNTimeBins->Multiply(fNLocalMaxima);
794     fNPads->Multiply(fNLocalMaxima);
795     fTimePosition->Multiply(fNLocalMaxima);
796     fIsAnalysed = kFALSE;
797   }
798
799   //
800   // TimeBin cut
801   //
802   if (iTimeBin<fFirstTimeBin) return 0;
803   if (iTimeBin>fLastTimeBin) return 0;
804   
805   // if pedestal calibrations are loaded subtract pedestals
806   if(fPedestal) {
807
808     Float_t ped = fPedestal->GetCalROC(iSector)->GetValue(iRow, iPad);
809     // Don't use data from pads where pedestals are abnormally small or large
810     if(ped<10 || ped>90)
811       return 0;
812     signal -= ped;
813   }
814   
815   if(fIsDQM) {
816
817     fOccVec->GetArray()[iSector] += 1.0;
818     // To change before committing
819     if(iPatch>=0 && iBranch>=0 && iPatch<=5 && iBranch <= 1)
820       fOccVecFine->GetArray()[(iSector%36)*12+iPatch*2+iBranch] += 1.0;
821   } else {
822     // In fNoThreshold we fill all data to estimate the ZS volume
823     Float_t count = fNoThreshold->GetCalROC(iSector)->GetValue(iRow, iPad);
824     fNoThreshold->GetCalROC(iSector)->SetValue(iRow, iPad,count+1);
825   }
826
827   // Require at least 3 ADC channels
828   if (signal < 3.0)
829     return 0;
830
831   // if noise calibrations are loaded require at least 3*sigmaNoise
832   if(fNoise) {
833     
834     Float_t noise = fNoise->GetCalROC(iSector)->GetValue(iRow, iPad);
835     
836     if(signal < noise*3.0)
837       return 0;
838   }
839
840   //
841   // This signal is ok and we store it in the cluster map
842   //
843
844   SetExpandDigit(iRow, iPad, iTimeBin, signal);
845
846   fSignalCounter++;
847   
848   return 1; // signal was accepted
849 }
850
851 //_____________________________________________________________________
852 void AliTPCdataQA::FindLocalMaxima(const Int_t iSector)
853 {
854   //
855   // This method is called after the data from each sector has been
856   // exapanded into an array
857   // Loop over the signals and identify local maxima and fill the
858   // calibration objects with the information
859   //
860
861   Int_t nLocalMaxima = 0;
862   const Int_t maxTimeBin = fTimeBinsMax+4; // Used to step between neighboring pads 
863                                            // Because we have tha pad-time data in a 
864                                            // 1d array
865
866   for (Int_t iRow = 0; iRow < fRowsMax; iRow++) {
867
868     Float_t* allBins = fAllBins[iRow];
869     Int_t* sigBins   = fAllSigBins[iRow];
870     const Int_t nSigBins   = fAllNSigBins[iRow];
871     
872     for (Int_t iSig = 0; iSig < nSigBins; iSig++) {
873
874       Int_t bin  = sigBins[iSig];
875       Float_t *b = &allBins[bin];
876
877       //
878       // Now we check if this is a local maximum
879       //
880
881       Float_t qMax = b[0];
882
883       // First check that the charge is bigger than the threshold
884       if (qMax<5) 
885         continue;
886       
887       // Require at least one neighboring pad with signal
888       if (b[-maxTimeBin]+b[maxTimeBin]<=0) continue;
889
890       // Require at least one neighboring time bin with signal
891       if (b[-1]+b[1]<=0) continue;
892       
893       //
894       // Check that this is a local maximum
895       // Note that the checking is done so that if 2 charges has the same
896       // qMax then only 1 cluster is generated 
897       // (that is why there is BOTH > and >=)
898       //
899       if (b[-maxTimeBin]   >= qMax) continue;
900       if (b[-1  ]          >= qMax) continue; 
901       if (b[+maxTimeBin]   > qMax)  continue; 
902       if (b[+1  ]          > qMax)  continue; 
903       if (b[-maxTimeBin-1] >= qMax) continue;
904       if (b[+maxTimeBin-1] >= qMax) continue; 
905       if (b[+maxTimeBin+1] > qMax)  continue; 
906       if (b[-maxTimeBin+1] >= qMax) continue;
907       
908       //
909       // Now we accept the local maximum and fill the calibration/data objects
910       //
911       nLocalMaxima++;
912
913       Int_t iPad, iTimeBin;
914       GetPadAndTimeBin(bin, iPad, iTimeBin);
915       
916       if(!fIsDQM) {
917         Float_t count = fNLocalMaxima->GetCalROC(iSector)->GetValue(iRow, iPad);
918         fNLocalMaxima->GetCalROC(iSector)->SetValue(iRow, iPad, count+1);
919         
920         count = fTimePosition->GetCalROC(iSector)->GetValue(iRow, iPad);
921         fTimePosition->GetCalROC(iSector)->SetValue(iRow, iPad, count+iTimeBin);
922       
923         Float_t charge = fMaxCharge->GetCalROC(iSector)->GetValue(iRow, iPad);
924         fMaxCharge->GetCalROC(iSector)->SetValue(iRow, iPad, charge + qMax);
925         
926         if(qMax>=10) {
927           count = fOverThreshold10->GetCalROC(iSector)->GetValue(iRow, iPad);
928           fOverThreshold10->GetCalROC(iSector)->SetValue(iRow, iPad, count+1);
929         }
930         if(qMax>=20) {
931           count = fOverThreshold20->GetCalROC(iSector)->GetValue(iRow, iPad);
932           fOverThreshold20->GetCalROC(iSector)->SetValue(iRow, iPad, count+1);
933         }
934         if(qMax>=30) {
935           count = fOverThreshold30->GetCalROC(iSector)->GetValue(iRow, iPad);
936           fOverThreshold30->GetCalROC(iSector)->SetValue(iRow, iPad, count+1);
937         }
938       }
939
940       //
941       // Calculate the total charge as the sum over the region:
942       //
943       //    o o o o o
944       //    o i i i o
945       //    o i C i o
946       //    o i i i o
947       //    o o o o o
948       //
949       // with qmax at the center C.
950       //
951       // The inner charge (i) we always add, but we only add the outer
952       // charge (o) if the neighboring inner bin (i) has a signal.
953       //
954       Int_t minP = 0, maxP = 0, minT = 0, maxT = 0;
955       Float_t qTot = qMax;
956       for(Int_t i = -1; i<=1; i++) {
957         for(Int_t j = -1; j<=1; j++) {
958           
959           if(i==0 && j==0)
960             continue;
961           
962           Float_t charge1 = GetQ(b, i, j, maxTimeBin, minT, maxT, minP, maxP);
963           qTot += charge1;
964           if(charge1>0) {
965             // see if the next neighbor is also above threshold
966             if(i*j==0) {
967               qTot += GetQ(b, 2*i, 2*j, maxTimeBin, minT, maxT, minP, maxP); 
968             } else {
969               // we are in a diagonal corner
970               qTot += GetQ(b,   i, 2*j, maxTimeBin, minT, maxT, minP, maxP); 
971               qTot += GetQ(b, 2*i,   j, maxTimeBin, minT, maxT, minP, maxP); 
972               qTot += GetQ(b, 2*i, 2*j, maxTimeBin, minT, maxT, minP, maxP); 
973             }
974           }
975         }
976       }
977       
978       if(fIsDQM) {
979         fHistQVsSector->Fill(iSector, qTot);
980         fHistQmaxVsSector->Fill(iSector, qMax);
981       } else {
982         Float_t charge = fMeanCharge->GetCalROC(iSector)->GetValue(iRow, iPad);
983         fMeanCharge->GetCalROC(iSector)->SetValue(iRow, iPad, charge + qTot);
984         
985         Float_t count = fNTimeBins->GetCalROC(iSector)->GetValue(iRow, iPad);
986         fNTimeBins->GetCalROC(iSector)->SetValue(iRow, iPad, count + maxT-minT+1);
987         
988         count = fNPads->GetCalROC(iSector)->GetValue(iRow, iPad);
989         fNPads->GetCalROC(iSector)->SetValue(iRow, iPad, count + maxP-minP+1);
990       
991         if((iSector%36)<18) { // A side
992           fHistQVsTimeSideA->Fill(iTimeBin, qTot);
993           fHistQMaxVsTimeSideA->Fill(iTimeBin, qMax);
994         } else {
995           fHistQVsTimeSideC->Fill(iTimeBin, qTot);
996           fHistQMaxVsTimeSideC->Fill(iTimeBin, qMax);      
997         }
998       }
999     } // end loop over signals
1000   } // end loop over rows
1001   
1002   fClusterCounter += nLocalMaxima;
1003 }
1004
1005 //_____________________________________________________________________
1006 void AliTPCdataQA::Analyse()
1007 {
1008   //
1009   //  Calculate calibration constants
1010   //
1011   
1012   AliInfo("Analyse called");
1013
1014   if(fIsDQM == kTRUE) {
1015     
1016     AliInfo("DQM flas is set -> No 2d information to analyze");
1017     return;
1018   }
1019
1020   if(fIsAnalysed == kTRUE) {
1021     
1022     AliInfo("No new data since Analyse was called last time");
1023     return;
1024   }
1025
1026   if(fEventCounter==0) {
1027     
1028       AliInfo("EventCounter == 0, Cannot analyse");
1029     return;
1030   }
1031   
1032   Int_t nTimeBins = fLastTimeBin - fFirstTimeBin +1;
1033   AliInfo(Form("EventCounter: %d , TimeBins: %d", fEventCounter, nTimeBins));
1034
1035   Float_t normalization = 1.0 / Float_t(fEventCounter * nTimeBins);
1036   fNoThreshold->Multiply(normalization);  
1037   
1038   fMeanCharge->Divide(fNLocalMaxima);
1039   fMaxCharge->Divide(fNLocalMaxima);
1040   fNTimeBins->Divide(fNLocalMaxima);
1041   fNPads->Divide(fNLocalMaxima);
1042   fTimePosition->Divide(fNLocalMaxima);
1043
1044   fIsAnalysed = kTRUE;
1045 }
1046
1047
1048 //_____________________________________________________________________
1049 void AliTPCdataQA::MakeTree(const char *fname) const {
1050   //
1051   // Export result to the tree -located in the file
1052   // This file can be analyzed using AliTPCCalibViewer
1053   // 
1054   AliTPCPreprocessorOnline preprocesor;
1055
1056   if (fNLocalMaxima) preprocesor.AddComponent(fNLocalMaxima);
1057   if (fMaxCharge) preprocesor.AddComponent(fMaxCharge);  
1058   if (fMeanCharge) preprocesor.AddComponent(fMeanCharge);  
1059   if (fNoThreshold) preprocesor.AddComponent(fNoThreshold);
1060   if (fNTimeBins) preprocesor.AddComponent(fNTimeBins);
1061   if (fNPads) preprocesor.AddComponent(fNPads);
1062   if (fTimePosition) preprocesor.AddComponent(fTimePosition);
1063   if (fOverThreshold10) preprocesor.AddComponent(fOverThreshold10);
1064   if (fOverThreshold20) preprocesor.AddComponent(fOverThreshold20);
1065   if (fOverThreshold30) preprocesor.AddComponent(fOverThreshold30);
1066
1067   preprocesor.DumpToFile(fname);  
1068 }
1069
1070
1071 //_____________________________________________________________________
1072 void AliTPCdataQA::MakeArrays(){
1073   //
1074   // The arrays for expanding the raw data are defined and 
1075   // som parameters are intialised
1076   //
1077   AliTPCROC * roc = AliTPCROC::Instance();
1078   //
1079   // To make the array big enough for all sectors we take 
1080   // the dimensions from the outer row of an OROC (the last sector)
1081   //
1082   fRowsMax     = roc->GetNRows(roc->GetNSector()-1);
1083   fPadsMax     = roc->GetNPads(roc->GetNSector()-1,fRowsMax-1);
1084   fTimeBinsMax = fLastTimeBin - fFirstTimeBin +1; 
1085
1086   //
1087   // Since we have added 2 pads (TimeBins) before and after the real pads (TimeBins) 
1088   // to make sure that we can always query the exanded table even when the 
1089   // max is on the edge
1090   //
1091
1092  
1093   fAllBins = new Float_t*[fRowsMax];
1094   fAllSigBins = new Int_t*[fRowsMax];
1095   fAllNSigBins = new Int_t[fRowsMax];
1096
1097   for (Int_t iRow = 0; iRow < fRowsMax; iRow++) {
1098     //
1099     Int_t maxBin = (fTimeBinsMax+4)*(fPadsMax+4);  
1100     fAllBins[iRow] = new Float_t[maxBin];
1101     memset(fAllBins[iRow],0,sizeof(Float_t)*maxBin); // set all values to 0
1102     fAllSigBins[iRow] = new Int_t[maxBin];
1103     fAllNSigBins[iRow] = 0;
1104   }
1105 }
1106
1107
1108 //_____________________________________________________________________
1109 void AliTPCdataQA::CleanArrays(){
1110   //
1111   //
1112   //
1113
1114   for (Int_t iRow = 0; iRow < fRowsMax; iRow++) {
1115
1116     // To speed up the performance by a factor 2 on cosmic data (and
1117     // presumably pp data as well) where the ocupancy is low, the
1118     // memset is only called if there is more than 1000 signals for a
1119     // row (of the order 1% occupancy)
1120     if(fAllNSigBins[iRow]<1000) {
1121       
1122       Float_t* allBins = fAllBins[iRow];
1123       Int_t* sigBins   = fAllSigBins[iRow];
1124       const Int_t nSignals = fAllNSigBins[iRow];
1125       for(Int_t i = 0; i < nSignals; i++)
1126         allBins[sigBins[i]]=0;      
1127     } else {
1128
1129       Int_t maxBin = (fTimeBinsMax+4)*(fPadsMax+4); 
1130       memset(fAllBins[iRow],0,sizeof(Float_t)*maxBin);
1131     }
1132
1133     fAllNSigBins[iRow]=0;
1134   }
1135 }
1136
1137 //_____________________________________________________________________
1138 void AliTPCdataQA::GetPadAndTimeBin(Int_t bin, Int_t& iPad, Int_t& iTimeBin){
1139   //
1140   // Return pad and timebin for a given bin
1141   //
1142   
1143   //  Int_t bin = iPad*(fTimeBinsMax+4)+iTimeBin;
1144   iTimeBin  = bin%(fTimeBinsMax+4);
1145   iPad      = (bin-iTimeBin)/(fTimeBinsMax+4);
1146
1147   iPad     -= 2;
1148   iTimeBin -= 2;
1149   iTimeBin += fFirstTimeBin;
1150
1151   R__ASSERT(iPad>=0 && iPad<=fPadsMax);
1152   R__ASSERT(iTimeBin>=fFirstTimeBin && iTimeBin<=fLastTimeBin);
1153 }
1154
1155 //_____________________________________________________________________
1156 void AliTPCdataQA::SetExpandDigit(const Int_t iRow, Int_t iPad, 
1157                                   Int_t iTimeBin, const Float_t signal) 
1158 {
1159   //
1160   // 
1161   //
1162   R__ASSERT(iRow>=0 && iRow<fRowsMax);
1163   R__ASSERT(iPad>=0 && iPad<=fPadsMax);
1164   R__ASSERT(iTimeBin>=fFirstTimeBin && iTimeBin<=fLastTimeBin);
1165   
1166   iTimeBin -= fFirstTimeBin;
1167   iPad     += 2;
1168   iTimeBin += 2;
1169   
1170   Int_t bin = iPad*(fTimeBinsMax+4)+iTimeBin;
1171   fAllBins[iRow][bin] = signal;
1172   fAllSigBins[iRow][fAllNSigBins[iRow]] = bin;
1173   fAllNSigBins[iRow]++;
1174 }
1175
1176 //______________________________________________________________________________
1177 Float_t AliTPCdataQA::GetQ(const Float_t* adcArray, const Int_t time, 
1178                            const Int_t pad, const Int_t maxTimeBins, 
1179                            Int_t& timeMin, Int_t& timeMax, 
1180                            Int_t& padMin,  Int_t& padMax) const
1181 {
1182   //
1183   // This methods return the charge in the bin time+pad*maxTimeBins
1184   // If the charge is above 0 it also updates the padMin, padMax, timeMin
1185   // and timeMax if necessary
1186   //
1187   Float_t charge = adcArray[time + pad*maxTimeBins];
1188   if(charge > 0) {
1189     timeMin = TMath::Min(time, timeMin); timeMax = TMath::Max(time, timeMax);
1190     padMin = TMath::Min(pad, padMin); padMax = TMath::Max(pad, padMax);
1191   }
1192   return charge; 
1193 }
1194
1195 //______________________________________________________________________________
1196 void AliTPCdataQA::Streamer(TBuffer &xRuub)
1197 {
1198   // Automatic schema evolution was first used from revision 4
1199   // Code based on:
1200   // http://root.cern.ch/root/roottalk/roottalk02/3207.html
1201
1202    UInt_t xRuus, xRuuc;
1203    if (xRuub.IsReading()) {
1204       Version_t xRuuv = xRuub.ReadVersion(&xRuus, &xRuuc);
1205       //we use the automatic algorithm for class version > 3
1206       if (xRuuv > 3) {
1207         AliTPCdataQA::Class()->ReadBuffer(xRuub, this, xRuuv, xRuus,
1208                                           xRuuc);
1209         return;
1210       }
1211       TH1F::Streamer(xRuub);
1212       xRuub >> fFirstTimeBin;
1213       xRuub >> fLastTimeBin;
1214       xRuub >> fAdcMin;
1215       xRuub >> fAdcMax;
1216       xRuub >> fNLocalMaxima;
1217       xRuub >> fMaxCharge;
1218       xRuub >> fMeanCharge;
1219       xRuub >> fNoThreshold;
1220       xRuub >> fNTimeBins;
1221       xRuub >> fNPads;
1222       xRuub >> fTimePosition;
1223       xRuub >> fEventCounter;
1224       xRuub >> fIsAnalysed;
1225       xRuub.CheckByteCount(xRuus, xRuuc, AliTPCdataQA::IsA());
1226    } else {
1227      AliTPCdataQA::Class()->WriteBuffer(xRuub,this);
1228    }
1229 }
1230
1231 //____________________________________________________________________________________________
1232 void AliTPCdataQA::FillOccupancyProfile()
1233 {
1234   // This has to be filled at the end of the loop over data
1235   if(!fIsDQM) 
1236     AliInfo("Method only meaningful for DQM");
1237   
1238   for(Int_t i = 0; i < 72; i++) {
1239
1240     fOccVec->GetArray()[i] /= fOccMaxVec->GetArray()[i];
1241     fHistOccVsSector->Fill(i, fOccVec->GetArray()[i]);
1242   }
1243
1244   const Int_t nBranches = 36*12;
1245   for(Int_t i = 0; i < nBranches; i++) {
1246
1247     fOccVecFine->GetArray()[i] /= fOccMaxVecFine->GetArray()[i];
1248
1249     const Int_t fullSector = Int_t(i/12);
1250
1251     Int_t branch = i - fullSector*12;
1252     const Int_t patch = Int_t(branch/2);
1253
1254     branch -= patch*2;
1255
1256     fHistOcc2dVsSector->Fill(fullSector+0.5*branch+0.1, patch+0.5, fOccVecFine->GetArray()[i]);
1257   }
1258 }
1259
1260 //____________________________________________________________________________________________
1261 void AliTPCdataQA::ResetProfiles()
1262 {
1263   if(!fIsDQM) 
1264     AliInfo("Method only meaningful for DQM");
1265   
1266   if(fHistQVsSector)
1267     fHistQVsSector->Reset();
1268   if(fHistQmaxVsSector)
1269     fHistQmaxVsSector->Reset();
1270   if(fHistOccVsSector)
1271     fHistOccVsSector->Reset();
1272   if(fHistOcc2dVsSector)
1273     fHistOcc2dVsSector->Reset();
1274
1275   if(fOccVec)
1276     for(Int_t i = 0; i < 72; i++)
1277       fOccVec->GetArray()[i] = 0.0;
1278   if(fOccVecFine)
1279     for(Int_t i = 0; i < 36*12; i++)
1280       fOccVecFine->GetArray()[i] = 0.0;
1281 }