]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCCalibRaw.cxx
style modifications (Markus)
[u/mrichter/AliRoot.git] / TPC / AliTPCCalibRaw.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 /* $Id:  $ */
17
18 ////////////////////////////////////////////////////////////////////////////////////////
19 //                                                                                    //
20 //     Implementation of the TPC Raw drift velocity and Altro L1 Phase  calibration   //
21 //                                                                                    //
22 //               Origin: Jens Wiechula, J.Wiechula@gsi.de                             //
23 //                                                                                    //
24 ////////////////////////////////////////////////////////////////////////////////////////
25 //
26 //
27 // *************************************************************************************
28 // *                                Class Description                                  *
29 // *************************************************************************************
30 /*
31
32 ----example---
33 TFile f("CalibAltro.root");
34 AliTPCCalibRaw *al=(AliTPCCalibRaw*)f.Get(f.GetListOfKeys()->At(0)->GetName())
35 {
36 TCanvas *c1=(TCanvas*)gROOT->FindObject("c1");
37 if (!c1) c1=new TCanvas("c1","c1");
38 c1->Clear();
39
40 TH2F h2f("h2","h2;RCU;fail",216,0,216,al->GetNevents(),0,al->GetNevents());
41 Bool_t first=kTRUE;
42 Int_t i,iev;
43 for (i=0;i<216;++i) {
44   TVectorF *v=al->GetALTROL1PhaseFailEventsRCU(i);
45   if (!v) continue;
46   for (iev=0;iev<al->GetNevents();++iev) {
47     h2f->SetBinContent(i+1,iev+1,(*v)(iev));
48   }
49 //   TH1F h(*v);
50 //   h.SetLineColor(i/216.*50+50);
51 //   ((TH1F*)h.Clone(Form("h%d",i)))->Draw(first?"":"same");
52 //   c1->Modified();
53 //   c1->Update();
54   first=kFALSE;
55 }
56 h2f->Draw("col");
57 }
58
59 */
60
61
62
63 //Root includes
64 #include <TROOT.h>
65 #include <TH2C.h>
66 #include <TH1F.h>
67 #include <TMap.h>
68 #include <TGraph.h>
69 #include <TObjString.h>
70 #include <TTimeStamp.h>
71 #include <TCanvas.h>
72 //AliRoot includes
73 #include "AliTPCCalROC.h"
74 #include "AliAltroRawStream.h"
75 #include "AliLog.h"
76 //class header
77 #include "AliTPCCalibRaw.h"
78
79 ClassImp(AliTPCCalibRaw)
80
81 AliTPCCalibRaw::AliTPCCalibRaw() :
82   AliTPCCalibRawBase(),
83   fPeakDetMinus(1),
84   fPeakDetPlus(2),
85   fNFailL1Phase(0),
86   fNFailL1PhaseEvent(0),
87   fFirstTimeStamp(0),
88   fNSecTime(600), //default 10 minutes
89   fNBinsTime(60), //default 60*10 minutes = 10 hours
90   fPadProcessed(kFALSE),
91   fCurrentChannel(-1),
92   fCurrentSector(-1),
93   fLastSector(-2),
94   fCurrentRow(-1),
95   fCurrentPad(-1),
96   fLastTimeBinProc(0),
97   fPeakTimeBin(0),
98   fLastSignal(0),
99   fNOkPlus(0),
100   fNOkMinus(0),
101   fNanoSec(0),
102   fArrCurrentPhaseDist(4),
103   fArrCurrentPhase(kNRCU),
104   fArrFailEventNumber(100),
105   fArrALTROL1Phase(100000),
106   fArrALTROL1PhaseEvent(kNRCU),
107   fArrALTROL1PhaseFailEvent(kNRCU),
108   fHnDrift(0x0),
109   fVOccupancyEvent(100000),
110   fVSignalSumEvent(100000),
111   fVOccupancySenEvent(100000),
112   fVSignalSumSenEvent(100000),
113   fVNfiredPadsSenEvent(100000),
114   fVTimeStampEvent(100000)
115 {
116   //
117   // Default ctor
118   //
119   SetNameTitle("AliTPCCalibRaw","AliTPCCalibRaw");
120   CreateDVhist();
121   for (Int_t ircu=0;ircu<kNRCU;++ircu) fArrCurrentPhase.GetMatrixArray()[ircu]=-1;
122   fFirstTimeBin=850;
123   fLastTimeBin=1020;
124 }
125 //_____________________________________________________________________
126 AliTPCCalibRaw::AliTPCCalibRaw(const TMap *config) :
127 AliTPCCalibRawBase(),
128 fPeakDetMinus(1),
129 fPeakDetPlus(2),
130 fNFailL1Phase(0),
131 fNFailL1PhaseEvent(0),
132 fFirstTimeStamp(0),
133 fNSecTime(600), //default 10 minutes
134 fNBinsTime(60), //default 60*10 minutes = 10 hours
135 fPadProcessed(kFALSE),
136 fCurrentChannel(-1),
137 fCurrentSector(-1),
138 fLastSector(-2),
139 fCurrentRow(-1),
140 fCurrentPad(-1),
141 fLastTimeBinProc(0),
142 fPeakTimeBin(0),
143 fLastSignal(0),
144 fNOkPlus(0),
145 fNOkMinus(0),
146 fNanoSec(0),
147 fArrCurrentPhaseDist(4),
148 fArrCurrentPhase(kNRCU),
149 fArrFailEventNumber(100),
150 fArrALTROL1Phase(100000),
151 fArrALTROL1PhaseEvent(kNRCU),
152 fArrALTROL1PhaseFailEvent(kNRCU),
153 fHnDrift(0x0),
154 fVOccupancyEvent(100000),
155 fVSignalSumEvent(100000),
156 fVOccupancySenEvent(100000),
157 fVSignalSumSenEvent(100000),
158 fVNfiredPadsSenEvent(100000),
159 fVTimeStampEvent(100000)
160 {
161   //
162   // Default ctor
163   //
164   SetNameTitle("AliTPCCalibRaw","AliTPCCalibRaw");
165   CreateDVhist();
166   for (Int_t ircu=0;ircu<kNRCU;++ircu) fArrCurrentPhase.GetMatrixArray()[ircu]=-1;
167   fFirstTimeBin=850;
168   fLastTimeBin=1020;
169   if (config->GetValue("FirstTimeBin")) fFirstTimeBin = ((TObjString*)config->GetValue("FirstTimeBin"))->GetString().Atoi();
170   if (config->GetValue("LastTimeBin")) fLastTimeBin = ((TObjString*)config->GetValue("LastTimeBin"))->GetString().Atoi();
171   if (config->GetValue("DebugLevel")) fDebugLevel = ((TObjString*)config->GetValue("DebugLevel"))->GetString().Atoi();
172 }
173
174 //_____________________________________________________________________
175 AliTPCCalibRaw::~AliTPCCalibRaw()
176 {
177   //
178   // dtor
179   //
180   delete fHnDrift;  
181 }
182 //_____________________________________________________________________
183 // AliTPCCalibRaw& AliTPCCalibRaw::operator = (const  AliTPCCalibRaw &source)
184 // {
185 //   //
186 //   // assignment operator
187 //   //
188 //   if (&source == this) return *this;
189 //   new (this) AliTPCCalibRaw(source);
190 //   
191 //   return *this;
192 // }
193
194 //_____________________________________________________________________
195 Int_t AliTPCCalibRaw::Update(const Int_t isector, const Int_t iRow, const Int_t iPad,
196              const Int_t iTimeBin, const Float_t signal)
197 {
198   //
199   // Data filling method
200   //
201   if (iRow<0) return 0;
202   if (iPad<0) return 0;
203   if (iTimeBin<0) return 0;
204   if (!fFirstTimeStamp) fFirstTimeStamp=GetTimeStamp();
205   //
206   Int_t iChannel  = fROC->GetRowIndexes(isector)[iRow]+iPad; //  global pad position in sector
207   //occupancy
208   fVOccupancyEvent.GetMatrixArray()[GetNevents()]++;
209   fVSignalSumEvent.GetMatrixArray()[GetNevents()]+=signal;
210   //occupancy in sensitive regions
211   Int_t npads=(Int_t)fROC->GetNPads(isector,iRow);
212   Int_t cpad=iPad-npads/2;
213   if (isector<(Int_t)fROC->GetNInnerSector()){
214     //IROC case (spot)
215     if ( iRow>19 && iRow<46 ){
216       if ( TMath::Abs(cpad)<7 ){
217         fVOccupancySenEvent.GetMatrixArray()[GetNevents()]++;
218         fVSignalSumSenEvent.GetMatrixArray()[GetNevents()]+=signal;
219         if (iChannel!=fCurrentChannel) fVNfiredPadsSenEvent.GetMatrixArray()[GetNevents()]++;
220       }
221     }
222   } else if ( iRow>75 ){
223     //OROC case (outer corners and last three rows are sensitive)
224     Int_t padEdge=(Int_t)TMath::Min(iPad,npads-iPad);
225     Int_t nrows=(Int_t)fROC->GetNRows(isector);
226     if ((nrows-iRow-1)<3 || padEdge<((((Int_t)iRow-76)/4+1))*2){
227       fVOccupancySenEvent.GetMatrixArray()[GetNevents()]++;
228       fVSignalSumSenEvent.GetMatrixArray()[GetNevents()]+=signal;
229       if (iChannel!=fCurrentChannel) fVNfiredPadsSenEvent.GetMatrixArray()[GetNevents()]++;
230     }
231   }
232   //
233   if ( (iTimeBin>fLastTimeBin) || (iTimeBin<fFirstTimeBin)   ) return 0;
234   //don't process edge pads
235   if (IsEdgePad(isector,iRow,iPad)) return 0;
236 //   Double_t x[kHnBinsDV]={1,isector,0};
237 //   fHnDrift->Fill(x);
238   if (fCurrentChannel==iChannel){
239     if (fPadProcessed) return 0;
240   } else {
241     fPadProcessed=kFALSE;
242     fNOkPlus=0;
243     fNOkMinus=0;
244     fPeakTimeBin=0;
245     fLastSignal=0;
246   }
247 //   Double_t x2[kHnBinsDV]={2,isector,0};
248 //   fHnDrift->Fill(x2);
249   
250
251   if (signal>fLastSignal) ++fNOkPlus;
252   else if(signal<fLastSignal && fNOkPlus>=fPeakDetPlus){
253     ++fNOkMinus;
254     if (!fPeakTimeBin) fPeakTimeBin=fLastTimeBinProc;
255     if ( fNOkMinus>=fPeakDetMinus ) {
256       Double_t x[kHnBinsDV]={fPeakTimeBin,isector,(fTimeStamp-fFirstTimeStamp)/fNSecTime};
257       fHnDrift->Fill(x);
258     } 
259   } else {
260     fNOkPlus=0;
261     fNOkMinus=0;
262     fPeakTimeBin=0;
263     fLastSignal=0;
264   }
265
266   fLastTimeBinProc=iTimeBin;
267   fLastSignal=TMath::Nint(signal);
268   fCurrentChannel = iChannel;
269   return 0;
270 }
271 //_____________________________________________________________________
272 void AliTPCCalibRaw::UpdateDDL(){
273   //
274   // fill ALTRO L1 information
275   //
276   
277   //set nanoseconds
278   if (!fNanoSec) {
279     TTimeStamp s;
280     fNanoSec=s.GetNanoSec();
281   }
282   // current phase
283   Int_t phase=(Int_t)(GetL1PhaseTB()*4.);
284   //Fill pahse information of current rcu and event
285   fArrCurrentPhase.GetMatrixArray()[fCurrDDLNum]=phase;
286   //increase phase counter
287   ++((fArrCurrentPhaseDist.GetMatrixArray())[phase]);
288   
289 }
290 //_____________________________________________________________________
291 void AliTPCCalibRaw::ResetEvent()
292 {
293   //
294   // Reset event counters
295   //
296
297   fCurrentChannel=-1;
298   fCurrentRow=-1;
299   fCurrentPad=-1;
300   fArrCurrentPhaseDist.Zero();
301 }
302 //_____________________________________________________________________
303 void AliTPCCalibRaw::EndEvent()
304 {
305   //
306   // End event analysis
307   //
308
309   
310   //find phase of the current event
311   Int_t phaseMaxEntries=-1;
312   Int_t maxEntries=0;
313   for (Int_t i=0;i<fArrCurrentPhaseDist.GetNrows();++i){
314     Int_t entries=(Int_t)fArrCurrentPhaseDist[i];
315     if (maxEntries<entries) {
316       maxEntries=entries;
317       phaseMaxEntries=i;
318     }
319   }
320   // store phase of current event
321   if (fArrALTROL1Phase.GetNrows()-1<=GetNevents())
322     fArrALTROL1Phase.ResizeTo(GetNevents()+10000);
323   (fArrALTROL1Phase.GetMatrixArray())[GetNevents()]=phaseMaxEntries;
324   
325   //loop over RCUs and test failures
326   UInt_t fail=0;
327   for (Int_t ircu=0;ircu<kNRCU;++ircu){
328     Int_t phase=(Int_t)fArrCurrentPhase[ircu];
329     if (phase<0) continue;
330     if (phase!=phaseMaxEntries){
331       TVectorF *arr=MakeArrL1PhaseRCU(fCurrDDLNum,kTRUE);
332       if (arr->GetNrows()-1<=(Int_t)fNFailL1PhaseEvent) arr->ResizeTo(arr->GetNrows()+100);
333       (arr->GetMatrixArray())[fNFailL1PhaseEvent]=phase;
334       ++fNFailL1Phase;
335       fail=1;
336       }
337     //reset current phase information
338     fArrCurrentPhase[ircu]=-1;
339   }
340   if (fail){
341     if (fArrFailEventNumber.GetNrows()-1<=(Int_t)fNFailL1PhaseEvent) fArrFailEventNumber.ResizeTo(fArrFailEventNumber.GetNrows()+100);
342     fArrFailEventNumber.GetMatrixArray()[fNFailL1PhaseEvent]=GetNevents();
343   }
344   fNFailL1PhaseEvent+=fail;
345   //time stamps
346   fVTimeStampEvent.GetMatrixArray()[GetNevents()]=GetTimeStamp()-fFirstTimeStamp+fNanoSec*1e-9;
347   fNanoSec=0;
348   //occupance related
349   if (fVOccupancyEvent.GetNrows()-1<=GetNevents()){
350     fVOccupancyEvent.ResizeTo(GetNevents()+10000);
351     fVSignalSumEvent.ResizeTo(GetNevents()+10000);
352     fVOccupancySenEvent.ResizeTo(GetNevents()+10000);
353     fVSignalSumSenEvent.ResizeTo(GetNevents()+10000);
354     fVTimeStampEvent.ResizeTo(GetNevents()+10000);
355     fVNfiredPadsSenEvent.ResizeTo(GetNevents()+10000);
356   }
357   IncrementNevents();
358 }
359 //_____________________________________________________________________
360 TH2C *AliTPCCalibRaw::MakeHistL1RCUEvents(Int_t type)
361 {
362   // Create a 2D histo RCU:Events indicating the there was a deviation
363   // from the mean L1 phase of the event
364   //
365   //type: 0=Failures, 1=Phases
366
367   //number of relavant events, depending on version
368   Int_t nevents=GetNevents();
369   //check version
370   Bool_t newVersion=kFALSE;
371   for (Int_t ircu=0; ircu<kNRCU; ++ircu){
372     const TVectorF *v=GetALTROL1PhaseEventsRCU(ircu);
373     if (!v) continue;
374     if ((UInt_t)(v->GetNrows())==fNFailL1PhaseEvent){
375       newVersion=kTRUE;
376       nevents=fNFailL1PhaseEvent;
377     }
378     break;
379   }
380   TH2C *h2 = new TH2C("hL1FailRCUEvents","L1 Failures;RCU;Event",kNRCU,0,kNRCU,nevents,0,nevents);
381   Int_t add=0;
382   for (Int_t ircu=0;ircu<kNRCU;++ircu) {
383     const TVectorF *v=GetALTROL1PhaseEventsRCU(ircu);
384     if (type==0){
385       add=1;
386       h2->SetMinimum(0);
387       h2->SetMaximum(2);
388     } else if (type==1) {
389       add=0;
390       h2->SetMinimum(0);
391       h2->SetMaximum(4);
392     }
393     if (!v) continue;
394     for (Int_t iev=0;iev<nevents;++iev) {
395       Float_t val=(*v)(iev);
396       Float_t phase=fArrALTROL1Phase.GetMatrixArray()[iev];
397       if (newVersion) {
398         Int_t event=(Int_t)fArrFailEventNumber.GetMatrixArray()[iev];
399         phase=fArrALTROL1Phase.GetMatrixArray()[event];
400       }
401       if (type==0) val=(val!=phase);
402       h2->SetBinContent(ircu+1,iev+1,val+add);
403     }
404   }
405   return h2;
406 }
407 //_____________________________________________________________________
408 TH1F *AliTPCCalibRaw::MakeHistL1PhaseDist()
409 {
410   //
411   // L1 phase distribution. Should be flat in ideal case
412   //
413   TH1F *h=new TH1F("L1phaseDist","Normalized L1 phase distribution;phase;fraction of events",4,0,4);
414   h->Sumw2();
415   for (Int_t iev=0;iev<GetNevents();++iev) h->Fill(fArrALTROL1Phase.GetMatrixArray()[iev]);
416   if (GetNevents()>0) h->Scale(1./GetNevents());
417   h->SetMinimum(0);
418   h->SetMaximum(1);
419   return h;
420 }
421 //_____________________________________________________________________
422 TVectorF *AliTPCCalibRaw::MakeVectL1PhaseDist()
423 {
424   //
425   // L1 phase distribution. Should be flat in ideal case
426   //
427   TVectorF *v=new TVectorF(4);
428   for (Int_t iev=0;iev<GetNevents();++iev) {
429     Int_t phase=(Int_t)fArrALTROL1Phase.GetMatrixArray()[iev];
430     ((v->GetMatrixArray())[phase])+=1./GetNevents();
431   }
432   return v;
433 }
434 //_____________________________________________________________________
435 TH2C *AliTPCCalibRaw::MakeHistL1RCUEventsIROC(Int_t type)
436 {
437   //
438   // Create a 2D histo RCU:Events indicating the there was a deviation
439   // from the mean L1 phase of the event
440   //
441   TH2C *h2 = new TH2C("hL1FailRCUEventsIROC","L1 Failures IROCs;RCU;Event",72,0,36,GetNevents(),0,GetNevents());
442   for (Int_t ircu=0;ircu<72;++ircu) {
443     const TVectorF *v=0;
444     if (type==0)      v=GetALTROL1PhaseFailEventsRCU(ircu);
445     else if (type==1) v=GetALTROL1PhaseEventsRCU(ircu);
446     if (!v) continue;
447     for (Int_t iev=0;iev<GetNevents();++iev) {
448       h2->SetBinContent(ircu+1,iev+1,(*v)(iev));
449     }
450   }
451   return h2;
452 }
453 //_____________________________________________________________________
454 TH2C *AliTPCCalibRaw::MakeHistL1RCUEventsOROC(Int_t type)
455 {
456   //
457   // Create a 2D histo RCU:Events indicating the there was a deviation
458   // from the mean L1 phase of the event
459   //
460   TH2C *h2 = new TH2C("hL1FailRCUEventsOROC","L1 Failures OROCs;RCU;Event",144,0,36,GetNevents(),0,GetNevents());
461   for (Int_t ircu=72;ircu<kNRCU;++ircu) {
462     const TVectorF *v=0;
463     if (type==0)      v=GetALTROL1PhaseFailEventsRCU(ircu);
464     else if (type==1) v=GetALTROL1PhaseEventsRCU(ircu);
465     if (!v) continue;
466     for (Int_t iev=0;iev<GetNevents();++iev) {
467       h2->SetBinContent(ircu-72+1,iev+1,(*v)(iev));
468     }
469   }
470   return h2;
471 }
472 //_____________________________________________________________________
473 void AliTPCCalibRaw::CreateDVhist()
474 {
475   //
476   // Setup the HnSparse for the drift velocity determination
477   //
478   if (fHnDrift) return;
479   //HnSparse bins
480   //time bin, roc, time
481   Int_t    bins[kHnBinsDV] = {fLastTimeBin-fFirstTimeBin, 72, fNBinsTime};
482   Double_t xmin[kHnBinsDV] = {fFirstTimeBin,0,0};
483   Double_t xmax[kHnBinsDV] = {fLastTimeBin,72,fNBinsTime};
484   fHnDrift=new THnSparseI("fHnDrift",Form("Drift velocity using last time bin;time bin[#times 100ns];ROC;Time bin [#times %us]",fNSecTime),kHnBinsDV, bins, xmin, xmax);
485     
486 }
487 //_____________________________________________________________________
488 void AliTPCCalibRaw::Analyse()
489 {
490   //
491   // Analyse Data
492   //
493
494   //resize arrays
495   fArrALTROL1Phase.ResizeTo(GetNevents());
496   for (Int_t ircu=0;ircu<kNRCU;++ircu){
497     TVectorF *arr=MakeArrL1PhaseRCU(ircu);//MakeArrL1PhaseRCU(ircu);
498     if (!arr) continue;
499     arr->ResizeTo(fNFailL1PhaseEvent);
500     fArrFailEventNumber.ResizeTo(fNFailL1PhaseEvent);
501 //    TVectorF *arrF=MakeArrL1PhaseFailRCU(ircu);
502 //     arrF->ResizeTo(1);
503   }
504   //resize occupancy arrays only save event occupancy in sensitive regions by default
505   //save the rest in debub mode
506   fVOccupancySenEvent.ResizeTo(GetNevents());
507   if (fDebugLevel>0){
508     fVOccupancyEvent.ResizeTo(GetNevents());
509     fVSignalSumEvent.ResizeTo(GetNevents());
510     fVSignalSumSenEvent.ResizeTo(GetNevents());
511     fVNfiredPadsSenEvent.ResizeTo(GetNevents());
512     fVTimeStampEvent.ResizeTo(GetNevents());
513   } else {
514     fVOccupancyEvent.ResizeTo(0);
515     fVSignalSumEvent.ResizeTo(0);
516     fVSignalSumSenEvent.ResizeTo(0);
517     fVNfiredPadsSenEvent.ResizeTo(0);
518     fVTimeStampEvent.ResizeTo(0);
519   }
520   //Analyse drift velocity TODO
521   
522 }
523 //_____________________________________________________________________
524 TGraph* AliTPCCalibRaw::MakeGraphOccupancy(const Int_t type, const Int_t xType)
525 {
526   //
527   // create occupancy graph (of samples abouve threshold)
528   // type=0:   number of samples
529   // type=1:   mean data volume (ADC counts/sample)
530   // type=2:   data volume (ADC counts)
531   // type=3:   samples per ADC count
532   // type=4:   sample occupancy
533   //
534   // type=5: number of sample sensitive / number of samples
535   //
536   // same in sensitive regions:
537   // type=10:  number of samples
538   // type=11:  mean data volume (ADC counts/sample)
539   // type=12:  data volume (ADC counts)
540   // type=13:  samples per ADC count
541   // type=14:  sample occupancy
542   //
543   // type=16: number of samples sensitive / number of pads sensitive
544   // type=17: pad occupancy in sensitive regions
545   // xType=0:  vs. time stamp
546   // xType=1:  vs. event counter
547   //
548
549   TString title("Event occupancy");
550   TString xTitle("Time");
551   TString yTitle("number of samples");
552   TGraph *gr=new TGraph(GetNevents());
553   if (fVSignalSumEvent.GetNrows()==0&&!(type==10||type==14)) return 0;
554   TVectorF *vOcc=&fVOccupancyEvent;
555   TVectorF *vSum=&fVSignalSumEvent;
556   TVectorF *vPads=&fVNfiredPadsSenEvent;
557   Double_t norm=557568.;
558   if (type!=14&&fVOccupancyEvent.GetNrows()==0){
559     AliWarning("In non debug mode only occupancy in sensitive regions vs. event awailable!!!");
560     return 0;
561   }
562   if (type>=10){
563     vOcc=&fVOccupancySenEvent;
564     vSum=&fVSignalSumSenEvent;
565     vPads=&fVNfiredPadsSenEvent;
566     norm=33012.;
567   }
568   for (Int_t i=0;i<GetNevents(); ++i){
569     Double_t nAboveThreshold=vOcc->GetMatrixArray()[i];
570     
571     Double_t nSumADC        =1;
572     Double_t timestamp      =1;
573     Double_t nPads          =1;
574
575     if (fVOccupancyEvent.GetNrows()>0){
576       nSumADC        =vSum->GetMatrixArray()[i];
577       timestamp      =fVTimeStampEvent.GetMatrixArray()[i]+fFirstTimeStamp;
578       nPads          =vPads->GetMatrixArray()[i];
579     }
580     Double_t x=timestamp;
581     Double_t y=0;
582     //
583     if (xType==1)     x=i;
584     //
585     if (type%10==0)            y=nAboveThreshold;
586     if (type%10==1&&nAboveThreshold>0)      y=nSumADC/nAboveThreshold;
587     if (type%10==2)            y=nSumADC;
588     if (type%10==3&&nSumADC>0) y=nAboveThreshold/nSumADC;
589     if (type%10==4)            y=nAboveThreshold/(norm*(fLastTimeBin-fFirstTimeBin));
590     if (type==5)               y=fVOccupancySenEvent.GetMatrixArray()[i]/fVOccupancyEvent.GetMatrixArray()[i];
591     if (type==16&&nPads>0)     y=nAboveThreshold/nPads;
592     if (type==17)              y=nPads/norm;
593     //
594     gr->SetPoint(i,x,y);
595   }
596   if (xType==1) xTitle="Event";
597   if (type%10==1)      yTitle="Mean ADC counts/sample";
598   else if (type%10==2) yTitle="Data volume [ADC counts]";
599   else if (type%10==3) yTitle="samples per ADC count";
600   else if (type%10==4) yTitle="sample occupancy";
601   if (type==5)         yTitle="N samples (sensitive) / N samples";
602   if (type%10==6)      yTitle="N samples / N pads";
603   if (type==17)        yTitle="Pad Occupancy";
604   if (type>=10)        yTitle+=" (sensitive)";
605   title=yTitle+":"+xTitle;
606   title+=";"+xTitle+";"+yTitle;
607   gr->SetTitle(title.Data());
608   gr->SetEditable(kFALSE);
609   return gr;
610 }
611 //_____________________________________________________________________
612 // TGraph* AliTPCCalibRaw::MakeGraphNoiseEvents()
613 // {
614   //
615   // Not implemented for the moment
616   //
617 //   return 0;  
618 // }
619 //_____________________________________________________________________
620 TCanvas* AliTPCCalibRaw::MakeCanvasOccupancy(const Int_t xType, Bool_t sen)
621 {
622   //
623   // Create a canvas with occupancy information of all 'type's (see MakeGraphOccupancy)
624   // xType=0: vs. timestamp
625   // xType=1: vs. event number
626   //
627   // sen=kTRUE: for sensitive regions
628   //
629
630   TString name("RawOccupancy_");
631   TString title("Raw Occupancy vs. ");
632   if (xType==0){
633     name+="Time";
634     title+="time";
635   } else if (xType==1){
636     name+="Event";
637     title+="event";
638   }
639   if (sen){
640     name+="Sen";
641     title+=" (sensitive)";
642   }
643   TCanvas *c=(TCanvas*)gROOT->GetListOfCanvases()->FindObject(name.Data());
644   if (!c) c=new TCanvas(name.Data(),title.Data());
645   c->Clear();
646   c->Divide(2,2);
647   for (Int_t i=0;i<4;++i){
648     c->cd(i+1);
649     TGraph *gr=MakeGraphOccupancy(i+10*(Int_t)sen,xType);
650     gr->Draw("alp");
651   }
652   return c;
653 }
654
655 //_____________________________________________________________________
656 void AliTPCCalibRaw::Merge(AliTPCCalibRaw * const sig)
657 {
658   //
659   // Merge sig with this instance
660   //
661
662   if (!sig) return;
663
664   //Add last time bin distribution histogram
665   fHnDrift->Add(sig->fHnDrift);
666
667   //Add occupancy data
668   
669 }
670
671 //_____________________________________________________________________
672 Long64_t AliTPCCalibRaw::Merge(TCollection * const list)
673 {
674   //
675   // Merge all objects of this type in list
676   //
677   
678   Long64_t nmerged=1;
679   
680   TIter next(list);
681   AliTPCCalibRaw *ce=0;
682   TObject *o=0;
683   
684   while ( (o=next()) ){
685     ce=dynamic_cast<AliTPCCalibRaw*>(o);
686     if (ce){
687       Merge(ce);
688       ++nmerged;
689     }
690   }
691   
692   return nmerged;
693 }
694