]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCMonitor.cxx
Adding alarms and graphical representation of alarms
[u/mrichter/AliRoot.git] / TPC / AliTPCMonitor.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 $Log$
18 Revision 1.3  2007/10/12 13:36:27  cvetan
19 Coding convention fixes from Stefan
20
21 Revision 1.2  2007/09/18 09:44:45  cvetan
22 Sorting out some issues concerning the compilation with and without DATE support
23
24 Revision 1.1  2007/09/17 10:23:31  cvetan
25 New TPC monitoring package from Stefan Kniege. The monitoring package can be started by running TPCMonitor.C macro located in macros folder.
26
27 */   
28
29 ////////////////////////////////////////////////////////////////////////
30 ////
31 //// AliTPCMonitor class
32 //// 
33 //// Main class for TPC Monitor
34 //// Monitor can handle rootified data, files and online streams in DATE format.
35 //// The monitor GUI is started by the macro TPCMonitor.C
36 //// 
37 //// In the The data are read in in two read cycles. 
38 //// If no sector is specified (sectorid==-1) all sectors are read and only the max value 
39 //// for each channel is stored in a global histogram for Side A and C.
40 //// In this way the whole TPC can be read in at once.
41 ////
42 //// If the sector is specified  only one sector is read in and additional quantities
43 //// e.g baseline and baseline rms are calculated and stored. 
44 //// 
45 //// Author: Stefan Kniege, IKF, Frankfurt
46 ////       
47 ////
48 /////////////////////////////////////////////////////////////////////////
49
50
51
52 #include "AliLog.h"
53 #include "AliTPCMonitor.h"   
54 #include "AliTPCMonitorMappingHandler.h"
55 #include "AliTPCMonitorFFT.h"
56 #include "AliRawReader.h"
57 #include "AliRawReaderRoot.h"
58 #include "AliRawEventHeaderBase.h"
59 #include "AliAltroRawStreamV3.h"
60 #include "TH2F.h" 
61 #include "TF1.h"
62 #include "TMath.h"
63 #include "TCanvas.h"
64 #include "TH3S.h"
65 #include "TLegend.h"
66 #include "TROOT.h" 
67 #include "TDirectory.h"
68 #include "TSystem.h"
69 #include "TPaveText.h"   
70 #include "TFile.h"
71 #include <Riostream.h>
72
73 ClassImp(AliTPCMonitor)
74
75 const Int_t AliTPCMonitor::fgkHwMaskFEC              = 0x0780;                          
76 const Int_t AliTPCMonitor::fgkHwMaskBranch           = 0x0800;                          
77 const Int_t AliTPCMonitor::fgkHwMaskFECChannel       = 0x007f;                          
78 const Int_t AliTPCMonitor::fgkHwMaskAltroChannel     = 0x000f;                          
79 const Int_t AliTPCMonitor::fgkHwMaskAltroChip        = 0x0070;                          
80 const Int_t AliTPCMonitor::fgkHwMaskRCU              = 0x7000;                          
81
82 //____________________________________________________________________________
83 AliTPCMonitor::AliTPCMonitor(char* name, char* title) : 
84 AliTPCMonitorConfig(name,title),
85 fPad(new Int_t*[GetMaxHwAddr()]),
86 fPadMapHw(new Float_t[GetMaxHwAddr()]),
87 fPadMapRCU(new Int_t*[GetMaxHwAddr()]),
88 fHistIROC(0),
89 fHistOROC(0),
90 fHistIROCIndex(0),
91 fHistOROCIndex(0),
92 fHistIROCTime(0),
93 fHistOROCTime(0),
94 fHistIROCClone(0),
95 fHistOROCClone(0),
96 fHistIROCRMS(0),
97 fHistOROCRMS(0),
98 fHistIROCBASE(0),
99 fHistOROCBASE(0),
100 fHistIROCSUM(0),
101 fHistOROCSUM(0),
102 fHistChannelTime(0),
103 fHistAddrMapIndex(0),
104 fHistAddrMaxAdc(0),
105 fHistAddrBaseMean(0),
106 fHistAddrMaxAdcX(0),
107 fHistAddrAdcSum(0),
108 fHistAddrBaseRms(0),
109 fHistDistrSumIROC(0),
110 fHistDistrMaxIROC(0),
111 fHistDistrSumOROC(0),
112 fHistDistrMaxOROC(0),
113 fHistDistrBase2dIROC(0),
114 fHistDistrBase2dOROC(0),
115 fHistDistrBaseRmsIROC(0),
116 fHistDistrBaseMeanIROC(0),
117 fHistDistrBaseRmsOROC(0),
118 fHistDistrBaseMeanOROC(0),
119 fHistGlobalMaxA(0),
120 fHistGlobalMaxC(0),
121 fHistList(new TObjArray()),
122
123 fkNRowsIroc(63),
124 fkNRowsOroc(96),
125
126 fkNPadsIroc(110),
127 fkNPadsOroc(140),
128 fkNPadMinIroc(-55),
129 fkNPadMinOroc(-70),
130 fkNPadMaxIroc(55),
131 fkNPadMaxOroc(70),
132 fVerb(0),
133 fLastEv(0),
134 fEventNumber(0),
135 fEventNumberOld(0),
136 fDisableFit(0),
137 fExecGlob(0),
138 fExecPlaneMax(0),
139 fExecPadIrocRms(0),
140 fExecPadOrocRms(0),
141 fRunId(0),
142 fEqId(0),
143 fPadUsedRoc(-1),
144 fPadUsedHwAddr(-1),
145 fGdcId(0),
146 fLdcId(0),
147 fLdcIdOld(1),
148 fMapEqidsSec(new Int_t*[36]),
149 fMapEqidsRcu(new Int_t[1000]),
150 fMirror(1),
151 fChannelIter(0),
152 fMapHand(0),
153 fRawReader(0)
154 {
155   // Constructor
156   
157   for(Int_t i = 0; i<GetMaxHwAddr(); i++) { fPad[i]       = new Int_t[GetTimeBins()];}
158   for(Int_t i = 0; i<GetMaxHwAddr(); i++) { fPadMapRCU[i] = new Int_t[6];}
159   for(Int_t i = 0; i<36; i++) { fMapEqidsSec[i]       = new Int_t[6];}
160   
161   if (gDirectory)
162   {
163     if (!gDirectory->GetList())
164     {
165       Warning("Build","Current directory is not a valid directory");
166       return;
167     }
168     AliTPCMonitor *hold = (AliTPCMonitor*)gDirectory->GetList()->FindObject(GetName());
169     if(hold)
170     {
171       Warning("Build","Replacing existing histogram: %s (Potential memory leak).",GetName());
172       gDirectory->GetList()->Remove(hold);
173     }
174     gDirectory->Append(this);
175   }
176   CreateHistos();
177   SetEqIds();
178   
179 }
180
181 //____________________________________________________________________________
182 AliTPCMonitor::AliTPCMonitor(const AliTPCMonitor &monitor):
183 AliTPCMonitorConfig(monitor.GetName(),monitor.GetTitle()),
184 fPad(new Int_t*[GetMaxHwAddr()]),
185 fPadMapHw(new Float_t[GetMaxHwAddr()]),
186 fPadMapRCU(new Int_t*[GetMaxHwAddr()]),
187 fHistIROC(0),
188 fHistOROC(0),
189 fHistIROCIndex(0),
190 fHistOROCIndex(0),
191 fHistIROCTime(0),
192 fHistOROCTime(0),
193 fHistIROCClone(0),
194 fHistOROCClone(0),
195 fHistIROCRMS(0),
196 fHistOROCRMS(0),
197 fHistIROCBASE(0),
198 fHistOROCBASE(0),
199 fHistIROCSUM(0),
200 fHistOROCSUM(0),
201 fHistChannelTime(0),
202 fHistAddrMapIndex(0),
203 fHistAddrMaxAdc(0),
204 fHistAddrBaseMean(0),
205 fHistAddrMaxAdcX(0),
206 fHistAddrAdcSum(0),
207 fHistAddrBaseRms(0),
208 fHistDistrSumIROC(0),
209 fHistDistrMaxIROC(0),
210 fHistDistrSumOROC(0),
211 fHistDistrMaxOROC(0),
212 fHistDistrBase2dIROC(0),
213 fHistDistrBase2dOROC(0),
214 fHistDistrBaseRmsIROC(0),
215 fHistDistrBaseMeanIROC(0),
216 fHistDistrBaseRmsOROC(0),
217 fHistDistrBaseMeanOROC(0),
218 fHistGlobalMaxA(0),
219 fHistGlobalMaxC(0),
220 fHistList(new TObjArray()),
221 fkNRowsIroc(monitor.fkNRowsIroc),
222 fkNRowsOroc(monitor.fkNRowsOroc),
223 fkNPadsIroc(monitor.fkNPadsIroc),
224 fkNPadsOroc(monitor.fkNPadsOroc),
225 fkNPadMinIroc(monitor.fkNPadMinIroc),
226 fkNPadMinOroc(monitor.fkNPadMinOroc),
227 fkNPadMaxIroc(monitor.fkNPadMaxIroc),
228 fkNPadMaxOroc(monitor.fkNPadMaxOroc),
229 fVerb(monitor.fVerb),
230 fLastEv(monitor.fLastEv),
231 fEventNumber(monitor.fEventNumber),
232 fEventNumberOld(monitor.fEventNumberOld),
233 fDisableFit(monitor.fDisableFit),
234 fExecGlob(monitor.fExecGlob),
235 fExecPlaneMax(monitor.fExecPlaneMax),
236 fExecPadIrocRms(monitor.fExecPadIrocRms),
237 fExecPadOrocRms(monitor.fExecPadOrocRms),
238 fRunId(monitor.fRunId),
239 fEqId(monitor.fEqId),
240 fPadUsedRoc(monitor.fPadUsedRoc),
241 fPadUsedHwAddr(monitor.fPadUsedHwAddr),
242 fGdcId(monitor.fGdcId),
243 fLdcId(monitor.fLdcId),
244 fLdcIdOld(monitor.fLdcIdOld),
245 fMapEqidsSec(new Int_t*[36]),
246 fMapEqidsRcu(new Int_t[1000]),
247 fMirror(monitor.fMirror),
248 fChannelIter(monitor.fChannelIter),
249 fMapHand(monitor.fMapHand),
250 fRawReader(monitor.fRawReader)
251 {
252   // copy constructor
253   
254   fHistIROC=(TH2F*)monitor.fHistIROC->Clone(); fHistList->Add(fHistIROC);
255   fHistOROC=(TH2F*)monitor.fHistOROC->Clone(); fHistList->Add(fHistOROC);
256   fHistIROCIndex=(TH2S*)monitor.fHistIROCIndex->Clone(); fHistList->Add(fHistIROCIndex);
257   fHistOROCIndex=(TH2S*)monitor.fHistOROCIndex->Clone(); fHistList->Add(fHistOROCIndex);
258   fHistIROCTime=(TH2F*)monitor.fHistIROCTime->Clone(); fHistList->Add(fHistIROCTime);
259   fHistOROCTime=(TH2F*)monitor.fHistOROCTime->Clone(); fHistList->Add(fHistOROCTime);
260   fHistIROCClone=(TH2F*)monitor.fHistIROCClone->Clone(); fHistList->Add(fHistIROCClone);
261   fHistOROCClone=(TH2F*)monitor.fHistOROCClone->Clone(); fHistList->Add(fHistOROCClone);
262   fHistIROCRMS=(TH2F*)monitor.fHistIROCRMS->Clone(); fHistList->Add(fHistIROCRMS);
263   fHistOROCRMS=(TH2F*)monitor.fHistOROCRMS->Clone(); fHistList->Add(fHistOROCRMS);
264   fHistIROCBASE=(TH2F*)monitor.fHistIROCBASE->Clone(); fHistList->Add(fHistIROCBASE);
265   fHistOROCBASE=(TH2F*)monitor.fHistOROCBASE->Clone(); fHistList->Add(fHistOROCBASE);
266   fHistIROCSUM=(TH2F*)monitor.fHistIROCSUM->Clone(); fHistList->Add(fHistIROCSUM);
267   fHistOROCSUM=(TH2F*)monitor.fHistOROCSUM->Clone(); fHistList->Add(fHistOROCSUM);
268   
269   fHistChannelTime=(TH2F*)monitor.fHistChannelTime->Clone(); fHistList->Add(fHistChannelTime);
270   
271   fHistAddrMapIndex=(TH1F*)monitor.fHistAddrMapIndex->Clone(); fHistList->Add(fHistAddrMapIndex);
272   fHistAddrMaxAdc=(TH1F*)monitor.fHistAddrMaxAdc->Clone(); fHistList->Add(fHistAddrMaxAdc);
273   fHistAddrBaseMean=(TH1F*)monitor.fHistAddrBaseMean->Clone(); fHistList->Add(fHistAddrBaseMean);
274   fHistAddrMaxAdcX=(TH1F*)monitor.fHistAddrMaxAdcX->Clone(); fHistList->Add(fHistAddrMaxAdcX);
275   fHistAddrAdcSum=(TH1F*)monitor.fHistAddrAdcSum->Clone(); fHistList->Add(fHistAddrAdcSum);
276   fHistAddrBaseRms=(TH1F*)monitor.fHistAddrBaseRms->Clone(); fHistList->Add(fHistAddrBaseRms);
277   
278   
279   fHistDistrSumIROC=(TH1F*)monitor.fHistDistrSumIROC->Clone(); fHistList->Add(fHistDistrSumIROC);
280   fHistDistrMaxIROC=(TH1F*)monitor.fHistDistrMaxIROC->Clone(); fHistList->Add(fHistDistrMaxIROC);
281   fHistDistrSumOROC=(TH1F*)monitor.fHistDistrSumOROC->Clone(); fHistList->Add(fHistDistrSumOROC);
282   fHistDistrMaxOROC=(TH1F*)monitor.fHistDistrMaxOROC->Clone(); fHistList->Add(fHistDistrMaxOROC);
283   
284   fHistDistrBase2dIROC=(TH2F*)monitor.fHistDistrBase2dIROC->Clone(); fHistList->Add(fHistDistrBase2dIROC);
285   fHistDistrBase2dOROC=(TH2F*)monitor.fHistDistrBase2dOROC->Clone(); fHistList->Add(fHistDistrBase2dOROC);
286   fHistDistrBaseRmsIROC=(TH1D*)monitor.fHistDistrBaseRmsIROC->Clone(); fHistList->Add(fHistDistrBaseRmsIROC);
287   fHistDistrBaseMeanIROC=(TH1D*)monitor.fHistDistrBaseMeanIROC->Clone(); fHistList->Add(fHistDistrBaseMeanIROC);
288   fHistDistrBaseRmsOROC=(TH1D*)monitor.fHistDistrBaseRmsOROC->Clone(); fHistList->Add(fHistDistrBaseRmsOROC);
289   fHistDistrBaseMeanOROC=(TH1D*)monitor.fHistDistrBaseMeanOROC->Clone(); fHistList->Add(fHistDistrBaseMeanOROC);
290   
291   fHistGlobalMaxA=(TH2S*)monitor.fHistGlobalMaxA->Clone(); fHistList->Add(fHistGlobalMaxA);
292   fHistGlobalMaxC=(TH2S*)monitor.fHistGlobalMaxC->Clone(); fHistList->Add(fHistGlobalMaxC);
293   
294   // fPad       = new Int_t*[monitor.GetMaxHwAddr()];
295   for(Int_t i = 0; i<GetMaxHwAddr(); i++) { fPad[i]       = new Int_t[monitor.GetTimeBins()];}
296   
297   //fPadMapRCU = new Int_t*[monitor.GetMaxHwAddr()];
298   for(Int_t i = 0; i<monitor.GetMaxHwAddr(); i++) { fPadMapRCU[i] = new Int_t[6];}
299   
300   //fPadMapHw  = new Float_t[monitor.GetMaxHwAddr()];
301   
302   //fMapEqidsRcu   = new Int_t[1000];
303   //fMapEqidsSec   = new Int_t*[36];
304   for(Int_t i = 0; i<36; i++) { fMapEqidsSec[i]       = new Int_t[6];}
305   SetEqIds();
306   
307   if (gDirectory)
308   {
309     if (!gDirectory->GetList())
310     {
311       Warning("Build","Current directory is not a valid directory");
312       
313     }
314     else
315     {
316       AliTPCMonitor *hold = (AliTPCMonitor*)gDirectory->GetList()->FindObject(monitor.GetName());
317       if(hold)
318       {
319         Warning("Build","Replacing existing histogram: %s (Potential memory leak).",monitor.GetName());
320         gDirectory->GetList()->Remove(hold);
321       }
322       gDirectory->Append(this);
323     }
324   }
325 }
326
327
328 //____________________________________________________________________________
329
330 AliTPCMonitor &AliTPCMonitor:: operator= (const AliTPCMonitor& monitor)
331 {
332   // assigment operator
333   if(this!=&monitor)
334   {
335     ((AliTPCMonitorConfig *)this)->operator=(monitor);
336     fkNRowsIroc=monitor.fkNRowsIroc;
337     fkNRowsOroc=monitor.fkNRowsOroc;
338     fkNPadsIroc=monitor.fkNPadsIroc;
339     fkNPadsOroc=monitor.fkNPadsOroc;
340     fkNPadMinIroc=monitor.fkNPadMinIroc;
341     fkNPadMinOroc=monitor.fkNPadMinOroc;
342     fkNPadMaxIroc=monitor.fkNPadMaxIroc;
343     fkNPadMaxOroc=monitor.fkNPadMaxOroc;
344     fVerb=monitor.fVerb;
345     fLastEv=monitor.fLastEv;
346     fEventNumber=monitor.fEventNumber;
347     fEventNumberOld=monitor.fEventNumberOld;
348     fDisableFit=monitor.fDisableFit;
349     fExecGlob=monitor.fExecGlob;
350     fExecPlaneMax=monitor.fExecPlaneMax;
351     fExecPadIrocRms=monitor.fExecPadIrocRms;
352     fExecPadOrocRms=monitor.fExecPadOrocRms;
353     fRunId=monitor.fRunId;
354     fEqId=monitor.fEqId;
355     fPadUsedRoc=monitor.fPadUsedRoc;
356     fPadUsedHwAddr=monitor.fPadUsedHwAddr;
357     fGdcId=monitor.fGdcId;
358     fLdcId=monitor.fLdcId;
359     fLdcIdOld=monitor.fLdcIdOld;
360     fMapHand=monitor.fMapHand;
361     fRawReader=monitor.fRawReader;
362     
363     
364     fHistList = new TObjArray();
365     fHistIROC=(TH2F*)monitor.fHistIROC->Clone(); fHistList->Add(fHistIROC);
366     fHistOROC=(TH2F*)monitor.fHistOROC->Clone(); fHistList->Add(fHistOROC);
367     fHistIROCIndex=(TH2S*)monitor.fHistIROCIndex->Clone(); fHistList->Add(fHistIROCIndex);
368     fHistOROCIndex=(TH2S*)monitor.fHistOROCIndex->Clone(); fHistList->Add(fHistOROCIndex);
369     fHistIROCTime=(TH2F*)monitor.fHistIROCTime->Clone(); fHistList->Add(fHistIROCTime);
370     fHistOROCTime=(TH2F*)monitor.fHistOROCTime->Clone(); fHistList->Add(fHistOROCTime);
371     fHistIROCClone=(TH2F*)monitor.fHistIROCClone->Clone(); fHistList->Add(fHistIROCClone);
372     fHistOROCClone=(TH2F*)monitor.fHistOROCClone->Clone(); fHistList->Add(fHistOROCClone);
373     fHistIROCRMS=(TH2F*)monitor.fHistIROCRMS->Clone(); fHistList->Add(fHistIROCRMS);
374     fHistOROCRMS=(TH2F*)monitor.fHistOROCRMS->Clone(); fHistList->Add(fHistOROCRMS);
375     fHistIROCBASE=(TH2F*)monitor.fHistIROCBASE->Clone(); fHistList->Add(fHistIROCBASE);
376     fHistOROCBASE=(TH2F*)monitor.fHistOROCBASE->Clone(); fHistList->Add(fHistOROCBASE);
377     fHistIROCSUM=(TH2F*)monitor.fHistIROCSUM->Clone(); fHistList->Add(fHistIROCSUM);
378     fHistOROCSUM=(TH2F*)monitor.fHistOROCSUM->Clone(); fHistList->Add(fHistOROCSUM);
379     
380     fHistChannelTime=(TH2F*)monitor.fHistChannelTime->Clone(); fHistList->Add(fHistChannelTime);
381     
382     fHistAddrMapIndex=(TH1F*)monitor.fHistAddrMapIndex->Clone(); fHistList->Add(fHistAddrMapIndex);
383     fHistAddrMaxAdc=(TH1F*)monitor.fHistAddrMaxAdc->Clone(); fHistList->Add(fHistAddrMaxAdc);
384     fHistAddrBaseMean=(TH1F*)monitor.fHistAddrBaseMean->Clone(); fHistList->Add(fHistAddrBaseMean);
385     fHistAddrMaxAdcX=(TH1F*)monitor.fHistAddrMaxAdcX->Clone(); fHistList->Add(fHistAddrMaxAdcX);
386     fHistAddrAdcSum=(TH1F*)monitor.fHistAddrAdcSum->Clone(); fHistList->Add(fHistAddrAdcSum);
387     fHistAddrBaseRms=(TH1F*)monitor.fHistAddrBaseRms->Clone(); fHistList->Add(fHistAddrBaseRms);
388     
389     
390     fHistDistrSumIROC=(TH1F*)monitor.fHistDistrSumIROC->Clone(); fHistList->Add(fHistDistrSumIROC);
391     fHistDistrMaxIROC=(TH1F*)monitor.fHistDistrMaxIROC->Clone(); fHistList->Add(fHistDistrMaxIROC);
392     fHistDistrSumOROC=(TH1F*)monitor.fHistDistrSumOROC->Clone(); fHistList->Add(fHistDistrSumOROC);
393     fHistDistrMaxOROC=(TH1F*)monitor.fHistDistrMaxOROC->Clone(); fHistList->Add(fHistDistrMaxOROC);
394     
395     fHistDistrBase2dIROC=(TH2F*)monitor.fHistDistrBase2dIROC->Clone(); fHistList->Add(fHistDistrBase2dIROC);
396     fHistDistrBase2dOROC=(TH2F*)monitor.fHistDistrBase2dOROC->Clone(); fHistList->Add(fHistDistrBase2dOROC);
397     fHistDistrBaseRmsIROC=(TH1D*)monitor.fHistDistrBaseRmsIROC->Clone(); fHistList->Add(fHistDistrBaseRmsIROC);
398     fHistDistrBaseMeanIROC=(TH1D*)monitor.fHistDistrBaseMeanIROC->Clone(); fHistList->Add(fHistDistrBaseMeanIROC);
399     fHistDistrBaseRmsOROC=(TH1D*)monitor.fHistDistrBaseRmsOROC->Clone(); fHistList->Add(fHistDistrBaseRmsOROC);
400     fHistDistrBaseMeanOROC=(TH1D*)monitor.fHistDistrBaseMeanOROC->Clone(); fHistList->Add(fHistDistrBaseMeanOROC);
401     
402     fHistGlobalMaxA=(TH2S*)monitor.fHistGlobalMaxA->Clone(); fHistList->Add(fHistGlobalMaxA);
403     fHistGlobalMaxC=(TH2S*)monitor.fHistGlobalMaxC->Clone(); fHistList->Add(fHistGlobalMaxC);
404     
405     fPad       = new Int_t*[monitor.GetMaxHwAddr()];
406     for(Int_t i = 0; i<GetMaxHwAddr(); i++) { fPad[i]       = new Int_t[monitor.GetTimeBins()];}
407     
408     fPadMapRCU = new Int_t*[monitor.GetMaxHwAddr()];
409     for(Int_t i = 0; i<monitor.GetMaxHwAddr(); i++) { fPadMapRCU[i] = new Int_t[6];}
410     
411     fPadMapHw  = new Float_t[monitor.GetMaxHwAddr()];
412     
413     fMapEqidsRcu   = new Int_t[1000];
414     fMapEqidsSec   = new Int_t*[36];
415     for(Int_t i = 0; i<36; i++) { fMapEqidsSec[i]       = new Int_t[6];}
416     SetEqIds();
417     
418     if (gDirectory)
419     {
420       if (!gDirectory->GetList())
421       {
422         Warning("Build","Current directory is not a valid directory");
423         return  *this;
424       }
425       AliTPCMonitor *hold = (AliTPCMonitor*)gDirectory->GetList()->FindObject(monitor.GetName());
426       if(hold)
427       {
428         Warning("Build","Replacing existing histogram: %s (Potential memory leak).",monitor.GetName());
429         gDirectory->GetList()->Remove(hold);
430       }
431       gDirectory->Append(this);
432     }
433   }
434   return *this;
435 }
436
437
438 //____________________________________________________________________________
439 AliTPCMonitor::~AliTPCMonitor() 
440 {
441   // Destructor
442   
443   for(Int_t i = 0; i<GetMaxHwAddr(); i++) { delete[] fPad[i] ;}
444   for(Int_t i = 0; i<GetMaxHwAddr(); i++) { delete[] fPadMapRCU[i];}
445   delete[] fPadMapHw ;
446   DeleteHistos();
447 }
448
449 //____________________________________________________________________________
450 void AliTPCMonitor::CreateHistos() 
451 {
452   // Create histograms to be displayed
453   
454   if(fVerb) cout << " create new ones " << endl;
455   fHistIROC            = new TH2F("fHistIROC"            ,"fHistIROC"           ,fkNRowsIroc,0,fkNRowsIroc,fkNPadsIroc, fkNPadMinIroc, fkNPadMaxIroc);       fHistList->Add(fHistIROC);
456   fHistOROC            = new TH2F("fHistOROC"            ,"fHistOROC"           ,fkNRowsOroc,0,fkNRowsOroc,fkNPadsOroc, fkNPadMinOroc, fkNPadMaxOroc);       fHistList->Add(fHistOROC);
457   
458   fHistIROCIndex       = new TH2S("fHistIROCIndex"       ,"fHistIROCIndex"      ,fkNRowsIroc,0,fkNRowsIroc,fkNPadsIroc, fkNPadMinIroc, fkNPadMaxIroc);       fHistList->Add(fHistIROCIndex);
459   fHistOROCIndex       = new TH2S("fHistOROCIndex"       ,"fHistOROCIndex"      ,fkNRowsOroc,0,fkNRowsOroc,fkNPadsOroc, fkNPadMinOroc, fkNPadMaxOroc);       fHistList->Add(fHistOROCIndex);
460   
461   fHistIROCTime        = new TH2F("fHistIROCTime"        ,"fHistIROCTime"       ,fkNRowsIroc,0,fkNRowsIroc,fkNPadsIroc, fkNPadMinIroc, fkNPadMaxIroc);       fHistList->Add(fHistIROCTime);
462   fHistOROCTime        = new TH2F("fHistOROCTime"        ,"fHistOROCTime"       ,fkNRowsOroc,0,fkNRowsOroc,fkNPadsOroc, fkNPadMinOroc, fkNPadMaxOroc);       fHistList->Add(fHistOROCTime);
463   
464   fHistIROCRMS         = new TH2F("fHistIROCRMS"         ,"fHistIROCRMS"        ,fkNRowsIroc,0,fkNRowsIroc,fkNPadsIroc, fkNPadMinIroc, fkNPadMaxIroc);       fHistList->Add(fHistIROCRMS);
465   fHistOROCRMS         = new TH2F("fHistOROCRMS"         ,"fHistOROCRMS"        ,fkNRowsOroc,0,fkNRowsOroc,fkNPadsOroc, fkNPadMinOroc, fkNPadMaxOroc);       fHistList->Add(fHistOROCRMS);
466   
467   fHistIROCSUM         = new TH2F("fHistIROCSUM"         ,"fHistIROCSUM"        ,fkNRowsIroc,0,fkNRowsIroc,fkNPadsIroc, fkNPadMinIroc, fkNPadMaxIroc);       fHistList->Add(fHistIROCSUM);
468   fHistOROCSUM         = new TH2F("fHistOROCSUM"         ,"fHistOROCSUM"        ,fkNRowsOroc,0,fkNRowsOroc,fkNPadsOroc, fkNPadMinOroc, fkNPadMaxOroc);       fHistList->Add(fHistOROCSUM);
469   
470   fHistIROCBASE        = new TH2F("fHistIROCBASE"        ,"fHistIROCBASE"       ,fkNRowsIroc,0,fkNRowsIroc,fkNPadsIroc, fkNPadMinIroc, fkNPadMaxIroc);       fHistList->Add(fHistIROCBASE);
471   fHistOROCBASE        = new TH2F("fHistOROCBASE"        ,"fHistOROCBASE"       ,fkNRowsOroc,0,fkNRowsOroc,fkNPadsOroc, fkNPadMinOroc, fkNPadMaxOroc);       fHistList->Add(fHistOROCBASE);
472   
473   
474   fHistChannelTime     = new TH2F("fHistChannelTime"     ,"fHistChannelTime"    ,GetNumOfChannels(),0,GetNumOfChannels(),GetTimeBins(),0,GetTimeBins());fHistList->Add(fHistChannelTime);
475   fHistAddrMapIndex    = new TH1F("fHistAddrMapIndex"    ,"fHistAddrMapIndex"   ,GetMaxHwAddr()    ,0,GetMaxHwAddr());                                  fHistList->Add(fHistAddrMapIndex);
476   fHistAddrMaxAdc      = new TH1F("fHistAddrMaxAdc"      ,"fHistAddrMaxAdc"     ,GetMaxHwAddr(),0,GetMaxHwAddr());                                      fHistList->Add(fHistAddrMaxAdc);
477   fHistAddrMaxAdcX     = new TH1F("fHistAddrMaxAdcX"     ,"fHistAddrMaxAdcX"    ,GetMaxHwAddr(),0,GetMaxHwAddr());                                      fHistList->Add(fHistAddrMaxAdcX);
478   fHistAddrBaseMean    = new TH1F("fHistAddrBaseMean"    ,"fHistAddrBaseMean"   ,GetMaxHwAddr(),0,GetMaxHwAddr());                                      fHistList->Add(fHistAddrBaseMean);
479   fHistAddrAdcSum      = new TH1F("fHistAddrAdcSum"      ,"fHistAddrAdcSum"     ,GetMaxHwAddr(),0,GetMaxHwAddr());                                      fHistList->Add(fHistAddrAdcSum);
480   fHistAddrBaseRms     = new TH1F("fHistAddrBaseRms"     ,"fHistAddrBaseRms"    ,GetMaxHwAddr(),0,GetMaxHwAddr());                                      fHistList->Add(fHistAddrBaseRms);
481   fHistDistrSumIROC    = new TH1F("fHistDistrSumIROC"    ,"fHistDistrSumIROC"   ,400,0.0,4000.0);                                                       fHistList->Add(fHistDistrSumIROC);
482   fHistDistrMaxIROC    = new TH1F("fHistDistrMaxIROC"    ,"fHistDistrMaxIROC"   ,500,0.0,1000.0);                                                       fHistList->Add(fHistDistrMaxIROC);
483   fHistDistrSumOROC    = new TH1F("fHistDistrSumOROC"    ,"fHistDistrSumOROC"   ,400,0.0,4000.0);                                                       fHistList->Add(fHistDistrSumOROC);
484   fHistDistrMaxOROC    = new TH1F("fHistDistrMaxOROC"    ,"fHistDistrMaxOROC"   ,500,0.0,1000.0);                                                       fHistList->Add(fHistDistrMaxOROC);
485   
486   fHistDistrBase2dIROC = new TH2F("fHistDistrBase2dIROC" ,"fHistDistrBase2dIROC",100,0.0,100.0,100,0.0,10.0);                                           fHistList->Add(fHistDistrBase2dIROC);
487   fHistDistrBase2dOROC = new TH2F("fHistDistrBase2dOROC" ,"fHistDistrBase2dOROC",100,0.0,100.0,100,0.0,10.0);                                           fHistList->Add(fHistDistrBase2dOROC);
488   
489   fHistGlobalMaxA      = new TH2S("SIDE A"               ,"SIDE A"              ,500,-3000,3000,500,-3000,3000);                                        fHistList->Add(fHistGlobalMaxA);
490   fHistGlobalMaxC      = new TH2S("SIDE C"               ,"SIDE C"              ,500,-3000,3000,500,-3000,3000);                                        fHistList->Add(fHistGlobalMaxC);
491   
492   ResetArrays();
493 }
494 //____________________________________________________________________________
495 Int_t AliTPCMonitor::ProcessEvent()
496 {
497   // Process Event
498   // Depending on the value of the sector id all sectors (sectorid == -1) are processed.
499   //
500   // In this case only the maximum values are calculated per pad and filled to the global histograms
501   // In a second loop the last processed(displayed) sector will be processed (sectorid!=-1)
502   // again and the baseline rms and further quantities are calculated
503   //
504   // If only one sector should be processed SetProcOneSector(1) should be set.
505   // In this case only the specified (last/last displayed) sector will be processed.
506   //
507   // If GetProcNextEvent()==0 the same event will be processed again
508   
509   
510   Int_t sectorid = 0;
511   Int_t retflag  = 0; // id of last sector + 1000,  or error flag
512   if(GetProcNextEvent()==1 && fLastEv) { AliInfo("Last event already processed"); }
513   if(GetProcNextEvent()==1) ResetSectorArray();
514   
515   
516   if(GetProcNextEvent()==0 || GetProcOneSector()==1 ) sectorid = GetLastSector();
517   else                                                sectorid = -1;
518   
519   // first iteration
520   retflag = ReadDataNew(sectorid);
521   
522   SetLastProcFile(GetFile());
523   
524   if(retflag>=10 && retflag<1000){ AliError("Could not read event properly: Check file name and format or try next event");  return 0  ;}
525   
526   DrawHists(3);
527   
528   // second iteration
529   if(sectorid==-1 && retflag >1000)
530   {
531     AliInfo("Second read cycle");
532     SetProcNextEvent(0);
533     if(GetLastSectorDisplayed()==-1) {sectorid =  GetLastSector()         ;    }
534     else                             {sectorid =  GetLastSectorDisplayed();  SetLastSector(sectorid)           ;  }
535     retflag = ReadDataNew(sectorid);
536   }
537   
538   SetLastSectorDisplayed(sectorid)  ;
539   fMapHand->ReadfecHwMap(GetLastSector());
540   FillHistsPadPlane();
541   DrawHists(1);
542   SetEventProcessed(1);
543   return retflag;
544 }
545 //__________________________________________________________________
546 Int_t AliTPCMonitor::ReadDataNew(Int_t secid)
547 {
548   // Read Data File/Stream  for specified Format.
549   // Payload will be extracted from either ROOT or DATE format
550   // and passed to FillHistsDecode for decoding of the adc information
551   
552   if (!fRawReader){
553     fRawReader = AliRawReader::Create(GetFile());
554     SetLastProcFile(GetFile());
555   } else if (strcmp(GetLastProcFile(),GetFile())!=0){
556     delete fRawReader;
557     fRawReader = AliRawReader::Create(GetFile());
558 //     printf("New file!!!\n");
559     SetLastProcFile(GetFile());
560   }
561   
562   if (!fRawReader){
563     AliWarning("Coult not initialize raw reader");
564     return 11;
565   }
566   
567   Bool_t skip=kTRUE;
568   while(skip && GetProcNextEvent())
569   {
570     if(fVerb) cout << "AliTPCMonitor::ReadDataNew get event " << endl;
571     if(fRawReader->IsA()==AliRawReaderRoot::Class()){
572       printf("Root, NextEvent: %d\n",GetEventID());
573       if (!fRawReader->GotoEvent(GetEventID())){AliError("Could not get next Event"); return 11 ;}
574     } else if(!fRawReader->NextEvent()) { AliError("Could not get next Event"); return 11 ;}
575     // skip all events but physics, calibration and software trigger events!
576     UInt_t eventType=fRawReader->GetType();
577     if ( !(eventType==AliRawEventHeaderBase::kPhysicsEvent ||
578            eventType==AliRawEventHeaderBase::kCalibrationEvent ||
579            eventType==AliRawEventHeaderBase::kSystemSoftwareTriggerEvent ||
580            eventType==AliRawEventHeaderBase::kDetectorSoftwareTriggerEvent) ) {
581              if (fVerb) cout<< "Skipping event! Its neither of 'physics, calibration and software trigger event'" << endl;
582              continue;
583            }
584     skip=kFALSE;
585     //test if the TPC has data
586     UChar_t *data=0;
587     fRawReader->Select("TPC");
588     if (!fRawReader->ReadNextData(data)) skip=kTRUE;
589     fEventNumber = fRawReader->GetEventIndex();
590     fEventNumberOld = fRawReader->GetEventIndex();
591   }
592   
593 //   printf("secid: %d\n",secid);
594   
595   //========================== Histogram filling ======================
596   ResetHistos() ;
597   AliAltroRawStreamV3 *altro=new AliAltroRawStreamV3(fRawReader);
598   altro->SelectRawData("TPC");
599   altro->Reset();
600   fChannelIter=0;
601   
602   Int_t     hw                = 0;
603   Int_t     nextHwAddress     = 0;
604   Int_t     rcupatch          = 0;
605   Int_t     maxADC            = 0;
606   Int_t     maxx              = 0;
607   Int_t     sum               = 0;
608   Int_t     sumn              = 0;
609   Int_t           lastrcuid   = 0;
610   
611   while ( altro->NextDDL() ){
612     fGdcId        = fRawReader->GetGDCId() ;
613     fLdcId        = fRawReader->GetLDCId() ;
614     fRunId        = fRawReader->GetRunNumber() ;
615     fEqId         = fRawReader->GetEquipmentId();
616     rcupatch      = GetRCUPatch(fRunId, fEqId);
617     Int_t rcupatchSector=rcupatch%6;
618     lastrcuid     = (rcupatch+1000);
619 //     printf("RCU patch: %d, LDC: %d, EqId: %d\n",rcupatch, fLdcId, fEqId);
620     if(fLdcIdOld!=fLdcId &&  fChannelIter!=0) {
621       if(secid==-1)
622       {
623         FillGlobal(GetLastSector());
624         ResetArrays();
625 //         printf("filled sector: %d\n",GetLastSector());
626         fChannelIter =0;
627       }
628       else
629       {
630         return lastrcuid;
631 //         if (rcupatch/6!=secid) continue;
632       }
633 //       fChannelIter=0;
634     }
635     if (!CheckEqId(secid,fEqId)) continue;
636     while ( altro->NextChannel() ){
637       hw=altro->GetHWAddress();
638       nextHwAddress         = ( hw + (rcupatchSector<<12) );
639       fPad[fChannelIter][0] = nextHwAddress ;
640       fPadMapHw[nextHwAddress]    =  fChannelIter  ;
641       maxADC=0;
642       while ( altro->NextBunch() ){
643         Int_t  startTbin    = (Int_t)altro->GetStartTimeBin();
644         Int_t  bunchlength  = (Int_t)altro->GetBunchLength();
645         const UShort_t *sig = altro->GetSignals();
646         for (Int_t iTimeBin = 0; iTimeBin<bunchlength; iTimeBin++){
647           Int_t adc=(Int_t)sig[iTimeBin];
648           Int_t ntime=startTbin-iTimeBin;
649           //fill channel information
650           fPad[fChannelIter][ntime]  = adc;
651           if( (adc>maxADC) && (ntime>=GetRangeMaxAdcMin()) && (ntime<GetRangeMaxAdcMax()  ))  {maxADC = adc;maxx = ntime ;}
652           if(              (ntime>=GetRangeSumMin())    && (ntime<GetRangeSumMax()     ))  {sum+=adc; sumn++;}
653           
654         }
655       }
656       //get pedestal, noise
657       Float_t pedestal=TMath::Mean(GetRangeBaseMax()-GetRangeBaseMin(),fPad[fChannelIter]+GetRangeBaseMin());
658       Float_t noise   =TMath::RMS(GetRangeBaseMax()-GetRangeBaseMin(),fPad[fChannelIter]+GetRangeBaseMin());
659       fHistAddrMaxAdc->SetBinContent(nextHwAddress,maxADC-GetPedestals()*pedestal);
660       
661       if(secid!=-1)
662       {
663         if(rcupatchSector<2)
664         {
665           fHistDistrBase2dIROC->Fill(pedestal,noise);
666           fHistDistrSumIROC->Fill(sum);
667           fHistDistrMaxIROC->Fill(maxADC-pedestal*GetPedestals());
668           fHistDistrSumIROC->Fill(sum -sumn*pedestal*GetPedestals());
669         }
670         else
671         {
672           fHistDistrBase2dOROC->Fill(pedestal,noise);
673           fHistDistrSumOROC->Fill(sum);
674           fHistDistrMaxOROC->Fill(maxADC-pedestal*GetPedestals());
675           fHistDistrSumOROC->Fill(sum -sumn*pedestal*GetPedestals());
676         }
677         
678         fHistAddrAdcSum->SetBinContent(  nextHwAddress,sum);
679         fHistAddrMapIndex->SetBinContent(nextHwAddress,fChannelIter);
680         fHistAddrBaseMean->SetBinContent(nextHwAddress,pedestal);
681         fHistAddrMaxAdcX->SetBinContent( nextHwAddress,maxx);
682         fHistAddrBaseRms->SetBinContent( nextHwAddress,noise);
683       }
684       
685       ++fChannelIter;
686     }
687     SetLastSector(rcupatch/6);
688     if(fChannelIter!=0 && secid==-1 ) SetSectorFilled(rcupatch/6);
689     fLdcIdOld = fLdcId ;
690   }
691   delete altro;
692   if(fChannelIter!=0 && secid==-1) { FillGlobal(GetLastSector());}
693   return lastrcuid;
694 }
695
696
697
698
699
700
701
702
703
704
705
706 //____________________________________________________________________________
707 void AliTPCMonitor::FillHistsPadPlane() 
708 {
709   // Fill 2Dim histograms for IROC and OROC (max , rms and sum)
710   
711   if(fVerb)cout << "AliTPCMonitor::FillHistsPadPlane() Start  " << endl;
712   if(fVerb)PrintConfig();
713   
714   Int_t pad    = 0;
715   Int_t row    = 0;
716   Int_t padmax = 0;
717   Int_t hwadd  = 0;
718   
719   for(Int_t ch = 0; ch<fChannelIter; ch++)
720   {
721     hwadd=  fPad[ch][0];
722     fHistChannelTime->SetCellContent(ch,0,hwadd);
723     
724     for(Int_t bin = 1; bin <GetTimeBins(); bin++)
725     {
726       if( fHistChannelTime->GetCellContent(ch,bin)!=0)  cout << " cellcontent already set " << endl;
727       if(     GetPedestals()==1 ) fHistChannelTime->SetCellContent(ch,bin,(fPad[ch][bin]- fHistAddrBaseMean->GetBinContent(hwadd)));
728       else                        fHistChannelTime->SetCellContent(ch,bin,(fPad[ch][bin]));
729     }
730     
731     pad    = fMapHand->GetPad(   hwadd);
732     row    = fMapHand->GetPadRow(hwadd);
733     padmax = fMapHand->GetNumofPads(row);
734     
735     if(row<63)
736     {
737       fHistIROC->SetCellContent(     row    +1 ,pad +55 -padmax/2 +1,fHistAddrMaxAdc->GetBinContent(  hwadd));
738       fHistIROCIndex->SetCellContent(row    +1 ,pad +55 -padmax/2 +1,ch);
739       fHistIROCRMS->SetCellContent(  row    +1 ,pad +55 -padmax/2 +1,fHistAddrBaseRms->GetBinContent( hwadd));
740       fHistIROCBASE->SetCellContent( row    +1 ,pad +55 -padmax/2 +1,fHistAddrBaseMean->GetBinContent(hwadd));
741       fHistIROCSUM->SetCellContent(  row    +1 ,pad +55 -padmax/2 +1,fHistAddrAdcSum->GetBinContent(  hwadd));
742     }
743     else
744     {
745       fHistOROC->SetCellContent(     row-63 +1 ,pad +70 -padmax/2 +1,fHistAddrMaxAdc->GetBinContent(  hwadd));
746       fHistOROCIndex->SetCellContent(row-63 +1 ,pad +70 -padmax/2 +1,ch);
747       fHistOROCRMS->SetCellContent(  row-63 +1 ,pad +70 -padmax/2 +1,fHistAddrBaseRms->GetBinContent( hwadd));
748       fHistOROCBASE->SetCellContent( row-63 +1 ,pad +70 -padmax/2 +1,fHistAddrBaseMean->GetBinContent(hwadd));
749       fHistOROCSUM->SetCellContent(  row-63 +1 ,pad +70 -padmax/2 +1,fHistAddrAdcSum->GetBinContent(  hwadd));
750     }
751   }
752   
753   fHistChannelTime->GetXaxis()->SetRange(0,fChannelIter);
754   fHistChannelTime->GetYaxis()->SetRange(0,GetTimeBins());
755 }
756
757
758
759 //____________________________________________________________________________
760 void AliTPCMonitor::ResetArrays() 
761 {
762   // Reset data arrays
763   for(Int_t row = 0 ; row < fkNRowsIroc; row++)
764   {
765     for(Int_t pad = 0 ; pad <  fkNPadsIroc ; pad++)
766     {
767       fHistIROCIndex->SetCellContent(row+1,pad+1,-1);
768     }
769   }
770   for(Int_t row = 0 ; row < fkNRowsOroc; row++)
771   {
772     for(Int_t pad = 0 ; pad <  fkNPadsOroc ; pad++)
773     {
774       fHistOROCIndex->SetCellContent(row+1,pad+1,-1);
775     }
776   }
777   
778   for(Int_t ch= 0; ch<GetMaxHwAddr(); ch++)
779   {
780     fHistAddrMaxAdcX->SetBinContent(ch,-1);
781     fHistAddrMapIndex->SetBinContent(ch,-1);
782     fHistAddrMaxAdc->SetBinContent(  ch, 0);
783     fHistAddrBaseMean->SetBinContent(  ch, 0);
784     fHistAddrAdcSum->SetBinContent(  ch, 0);
785     fHistAddrBaseRms->SetBinContent(ch, 0);
786   }
787   
788   for(Int_t ch = 0; ch< GetNumOfChannels(); ch++)
789   {
790     for(Int_t bin = 0; bin< GetTimeBins(); bin++)
791     {
792       fPad[ch][bin] = 0;
793     }
794   }
795   for(Int_t ch = 0; ch< GetMaxHwAddr(); ch++)
796   {
797     fPadMapHw[ch]=-1;
798     fPadMapRCU[ch][0]=-1;
799     fPadMapRCU[ch][1]=-1;
800     fPadMapRCU[ch][2]=-1;
801     fPadMapRCU[ch][3]=-1;
802     fPadMapRCU[ch][4]=-1;
803     fPadMapRCU[ch][5]=-1;
804   }
805   fPadUsedHwAddr=-1;
806 }
807
808
809 //____________________________________________________________________________
810 void AliTPCMonitor::ResetHistos() 
811 {
812   // Reset all but
813   for(Int_t i =0; i<fHistList->GetEntries(); i++)
814   {
815     if(GetProcNextEvent()==0 && strcmp(((TH1*)fHistList->At(i))->GetName(),"SIDE A")==0) continue;
816     if(GetProcNextEvent()==0 && strcmp(((TH1*)fHistList->At(i))->GetName(),"SIDE C")==0) continue;
817     ((TH1*)fHistList->At(i))->Reset();
818   }
819   ResetArrays();
820 }
821
822 //____________________________________________________________________________
823 void AliTPCMonitor::DeleteHistos() 
824 {
825   // Delete histograms
826   for(Int_t i =0; i<fHistList->GetEntries(); i++)
827   {
828     delete (TH1*)fHistList->At(i);
829   }
830 }
831
832
833 //__________________________________________________________________
834 Int_t AliTPCMonitor::CheckEqId(Int_t secid,Int_t eqid)
835 {
836   // Check if equipment id corresponds to any  rcu patch in sector
837   // Equipment ids changed during commisioning in 2006 (starting from run 704)
838   // However Runids started from 0 again in 2007
839   // Now only runids from commissioning in 2006 after runid 704 and all new ones are supported.
840   // Comment in equipment check for runids < 704 if old runs should be processed
841   
842   if(fVerb) cout << "AliTPCMonitor::CheckEqId  : SectorId  " << secid << " EquipmentId " << eqid  << " runid " << fRunId <<  endl;
843   Int_t passed =1;
844   //skip all eqids which do not belong to the TPC
845   if ( eqid<768||eqid>983 ) return 0;
846   //
847   if(fRunId<704 && 0) // commented out --> runs with runid < 704 in 2006 are not recognized anymore
848   {
849     if( (secid>-1) && (secid<36) )   // if ( secid is in range) { take only specific eqids}  else { take all }
850     {
851       if(      (secid==13) && ( eqid!=408 && eqid!=409  &&  eqid!=509 && eqid!=512  && eqid!=513 && eqid!=517 )) {passed=0;}
852       else if( (secid==4)  && ( eqid!=404 && eqid!=504  &&  eqid!=407 && eqid!=503  && eqid!=508 && eqid!=506 )) {passed=0;}
853     }
854     else                                                                   {if(fVerb) cout << "passed check "<< endl; }
855   }
856   else
857   {
858     if( (secid>-1) && (secid<36) )   // if ( secid is in range) { take only specific eqids}  else { take all }
859     {
860       if(eqid!=fMapEqidsSec[secid][0] && eqid!= fMapEqidsSec[secid][1] && eqid!=fMapEqidsSec[secid][2] &&
861          eqid!=fMapEqidsSec[secid][3] && eqid!= fMapEqidsSec[secid][4] && eqid!=fMapEqidsSec[secid][5] )  {passed=0;}
862     }
863     else                                                                   {if(fVerb) cout << "passed check "<< endl;}
864   }
865   
866   return passed;
867 }
868
869 //__________________________________________________________________
870 void AliTPCMonitor::SetEqIds()
871 {
872   // Set mapping for equipment ids
873   for(Int_t i = 0; i<36 ; i++)
874   {
875     for(Int_t j = 0; j<6; j++)
876     {
877       if(j<2) fMapEqidsSec[i][j]= 768+i*2+j;
878       else    fMapEqidsSec[i][j]= 840+i*4+j-2;
879     }
880   }
881   
882   for(Int_t i = 0; i<36 ; i++)
883   {
884     for(Int_t j = 0; j<6; j++)
885     {
886       if(j<2) fMapEqidsRcu[768+i*2+j]   = i*6 +j;
887       else    fMapEqidsRcu[840+i*4+j-2] = i*6 +j;
888     }
889   }
890 }
891
892 //__________________________________________________________________
893 void AliTPCMonitor::FillGlobal(Int_t sector) 
894 {
895   
896   // Fill global histograms with max adc for each channel
897   
898   TH2S* hglob =0;
899   if((sector/18) ==0) hglob =  fHistGlobalMaxA;
900   else                hglob =  fHistGlobalMaxC;
901   
902   Float_t  rotsec  = (2*TMath::Pi()/18.0);
903   Float_t  rot     = (-rotsec*(sector%18)  +4*rotsec);
904   
905   Float_t  m11     =    TMath::Cos(rot);
906   Float_t  m12     =    TMath::Sin(rot);
907   Float_t  m21     = -1*TMath::Sin(rot);
908   Float_t  m22     =    TMath::Cos(rot);
909   
910   Int_t    max     = 0; // use integer for global view
911   
912   Double_t xval    = 0.0;
913   Double_t yval    = 0.0;
914   
915   Int_t    pad     = 0;
916   Int_t    row     = 0;
917   Int_t    padmax  = 0;
918   
919   Float_t  xdr     = 0;
920   Float_t  ydr     = 0;
921   
922   for(Int_t hw = 0; hw<fHistAddrMaxAdc->GetNbinsX(); hw++)
923   {
924     max = (Int_t)fHistAddrMaxAdc->GetBinContent(hw);
925     if(max!=-1)
926     {
927       pad    = fMapHand->GetPad(   hw);
928       row    = fMapHand->GetPadRow(hw);
929       padmax = fMapHand->GetNumofPads(row);
930       if (sector%36>17) fMirror=-1;
931       else fMirror=1;
932       GetXY(xval ,yval , padmax,row ,pad);
933       xdr    =  xval*m11 +yval*m12;
934       ydr    =  xval*m21 +yval*m22;
935       if(hglob->GetBinContent(hglob->GetXaxis()->FindBin(xdr),hglob->GetYaxis()->FindBin(ydr))==0)   hglob->Fill(xdr,ydr,(Int_t)max);
936     }
937   }
938 }
939
940
941 //__________________________________________________________________
942 void AliTPCMonitor::GetXY( Double_t& xval , Double_t& yval , Int_t padmax, Int_t row , Int_t pad) const 
943 {
944   // Get x and y position of pad
945   
946   if(row<63)
947   {
948     xval = fMirror*( 2*padmax -4*pad -2);
949     yval = 852.25 +7.5*row;
950   }
951   else
952   {
953     xval = fMirror*( 3*padmax -6*pad -3);
954     if((row-63)<63) {     yval = 10*(row-63) +1351;     }
955     else              {   yval = 15*(row-63-64)+1993.5; }
956   }
957   
958
959
960 //__________________________________________________________________
961 Int_t AliTPCMonitor::GetPadAtX(Float_t xval, Int_t row, Int_t padmax) const
962 {
963   // Get pad number at given position in x
964   Int_t pad    = 0;
965   
966   if(row<63) {pad = (Int_t)( ( (xval/fMirror) +2 -2*padmax)/-4);}
967   else       {pad = (Int_t)( ( (xval/fMirror) +3 -3*padmax)/-6);}
968   
969   if(pad>=padmax) return -1;
970   else            return pad ;
971   
972 }
973
974 //__________________________________________________________________
975 Int_t AliTPCMonitor::GetPadAtX(Float_t xval, Int_t row) const
976 {
977   
978  // Get pad number at given position in x
979   
980   Int_t padmax = fMapHand->GetNumofPads(row);
981   Int_t pad    = 0;
982   
983   if(row<63) {pad = (Int_t)( ( (xval/fMirror) +2 -2*padmax)/-4);}
984   else       {pad = (Int_t)( ( (xval/fMirror) +3 -3*padmax)/-6);}
985   
986   if(pad>=padmax) return -1;
987   else            return pad ;
988   
989 }
990
991 //__________________________________________________________________
992 void AliTPCMonitor::DrawHists(Int_t histos) 
993 {
994   
995   // Draw sets of histograms
996   // histos==1 : 2Dim histos for MAX adc  and add executables
997   // histos==2 : distributions max/rms/sum
998   // histos==3 : global max adc for specified SideA/C
999   
1000   
1001   if(fVerb)    cout << " Draw histos " << endl;
1002   Char_t cside[10];
1003   if(GetLastSector()/18==0 ) sprintf(cside,"A");
1004   else                       sprintf(cside,"C");
1005   
1006   Char_t titleSEC[256];   sprintf(titleSEC   ,"Sector %i Side %s Run : %05i EventID %i "       ,GetLastSector()%18,cside,fRunId, fEventNumber);
1007   Char_t titleEvent[256]; sprintf(titleEvent ,"Time <-> Channles  %s"                          ,titleSEC);
1008   Char_t titleIROC[256];  sprintf(titleIROC  ,"IROC %s"                                        ,titleSEC);
1009   Char_t titleOROC[256];  sprintf(titleOROC  ,"OROC %s"                                        ,titleSEC);
1010   
1011   Char_t titleMAX[256];   sprintf(titleMAX   ,"Max (timebin: %i,%i) %s"                        ,GetRangeMaxAdcMin(),GetRangeMaxAdcMax(),titleSEC);
1012   Char_t titleSUM[256];   sprintf(titleSUM   ,"Sum (timebin: %i,%i) %s"                        ,GetRangeSumMin()   ,GetRangeSumMax()   ,titleSEC);
1013   Char_t titleBASE[256];  sprintf(titleBASE  ,"Baseline RMS<->Mean  (timebin: %i-%i) %s"       ,GetRangeBaseMin()  ,GetRangeBaseMax()  ,titleSEC);
1014   Char_t titleMEAN[256];  sprintf(titleMEAN  ,"Baseline Mean (timebin: %i-%i) %s"              ,GetRangeBaseMin()  ,GetRangeBaseMax()  ,titleSEC);
1015   Char_t titleRMS[256] ;  sprintf(titleRMS   ,"Baseline RMS (timebin: %i-%i) %s"               ,GetRangeBaseMin()  ,GetRangeBaseMax()  ,titleSEC);
1016   
1017   if(histos==1)
1018   {
1019       // IROC _______________________________________________________________
1020     TCanvas* ciroc = 0;
1021     ciroc = (TCanvas*)gROOT->GetListOfCanvases()->FindObject("ciroc");
1022     if(!ciroc)
1023     {
1024       ciroc = CreateCanvas("ciroc");
1025       fExecPlaneMax=0;
1026     }
1027     ciroc->cd();
1028     
1029     fHistIROC->SetXTitle("row");
1030     fHistIROC->SetYTitle("pad");
1031     if(GetPedestals()) fHistIROC->SetZTitle("max ADC (baseline sub)");
1032     else               fHistIROC->SetZTitle("max ADC ");
1033     fHistIROC->SetTitle(titleIROC);
1034     fHistIROC->SetMinimum(0.01);
1035     fHistIROC->Draw("COLZ");
1036     ciroc->UseCurrentStyle();
1037     
1038     
1039     fHistIROCTime->SetXTitle("row"); fHistIROCTime->SetZTitle("peak time (fit)");       fHistIROCTime->SetYTitle("pad");      fHistIROCTime->SetTitle(titleIROC);
1040     fHistIROCRMS->SetXTitle("row");  fHistIROCRMS->SetZTitle( "baseline rms (ADC)");    fHistIROCRMS->SetYTitle("pad");       fHistIROCRMS->SetTitle(titleIROC);
1041     
1042       // OROC
1043     TCanvas* coroc = 0;
1044     coroc =(TCanvas*)gROOT->GetListOfCanvases()->FindObject("coroc");
1045     if(!coroc) {
1046       coroc = CreateCanvas("coroc");
1047       fExecPlaneMax=0;
1048     }
1049     coroc->cd();
1050     
1051     fHistOROC->SetXTitle("row");
1052     fHistOROC->SetYTitle("pad");
1053     if(GetPedestals()) fHistOROC->SetZTitle("max ADC (baseline sub)");
1054     else               fHistOROC->SetZTitle("max ADC ");
1055     fHistOROC->SetTitle(titleOROC);
1056     fHistOROC->SetMinimum(0.01);
1057     fHistOROC->Draw("COLZ");
1058     coroc->UseCurrentStyle();
1059     
1060     
1061     fHistOROCTime->SetXTitle("row"); fHistOROCTime->SetZTitle("peak time (fit) (timebins)"); fHistOROCTime->SetYTitle("pad"); fHistOROCTime->SetTitle(titleOROC);
1062     fHistOROCRMS->SetXTitle("row");  fHistOROCRMS->SetZTitle("baseline rms (ADC)");          fHistOROCRMS->SetYTitle("pad");  fHistOROCRMS->SetTitle(titleOROC);
1063     
1064       // SUM
1065     Char_t namesum[256] ; sprintf(namesum,"ADC sum (bins: %i, %i)",GetRangeSumMin() ,GetRangeSumMax() );
1066     fHistIROCSUM->SetXTitle("row");      fHistIROCSUM->SetZTitle(namesum);      fHistIROCSUM->SetYTitle("pad");      fHistIROCSUM->SetTitle(titleIROC);
1067     fHistOROCSUM->SetXTitle("row");      fHistOROCSUM->SetZTitle(namesum);      fHistOROCSUM->SetYTitle("pad");      fHistOROCSUM->SetTitle(titleOROC);
1068     
1069       // BASE
1070     Char_t namebase[256] ; sprintf(namebase ,"base mean (timbebin: %i, %i )",GetRangeBaseMin(),GetRangeBaseMax());
1071     fHistIROCBASE->SetXTitle("row"); fHistIROCBASE->SetZTitle(namebase);  fHistIROCBASE->SetYTitle("pad");      fHistIROCBASE->SetTitle(titleIROC);
1072     fHistOROCBASE->SetXTitle("row"); fHistOROCBASE->SetZTitle(namebase);  fHistOROCBASE->SetYTitle("pad");      fHistOROCBASE->SetTitle(titleOROC);
1073     
1074     if(fHistIROCClone) fHistIROCClone->Delete();
1075     if(fHistOROCClone) fHistOROCClone->Delete();
1076     fHistIROCClone = (TH2F*)fHistIROC->Clone("fHistIROCClone");
1077     fHistOROCClone = (TH2F*)fHistOROC->Clone("fHistOROCClone");
1078     
1079       // Executables
1080     if(fExecPlaneMax==0)
1081     {
1082       Char_t carry1[100];
1083       sprintf(carry1,".x %s/TPC/AliTPCMonitorExec.C(1)",gSystem->Getenv("ALICE_ROOT"));
1084       ciroc->AddExec("pad",carry1);
1085       coroc->AddExec("pad",carry1);
1086       
1087       Char_t carry2[100];
1088       sprintf(carry2,".x %s/TPC/AliTPCMonitorExec.C(2)",gSystem->Getenv("ALICE_ROOT"));
1089       ciroc->AddExec("row",carry2);
1090       coroc->AddExec("row",carry2);
1091       fExecPlaneMax=1;
1092     }
1093     coroc->Update();
1094     ciroc->Update();
1095   }
1096   else if(histos==2)
1097   {
1098       // MAX ADC distribution  ____________________________________________
1099     TCanvas* cmax = 0;
1100     cmax = (TCanvas*)gROOT->GetListOfCanvases()->FindObject("cmax");
1101     if(!cmax)  cmax = CreateCanvas("cmax");
1102     
1103     cmax->cd();
1104     fHistDistrMaxIROC->GetXaxis()->SetRangeUser(0.0,1000.0);
1105     fHistDistrMaxIROC->SetXTitle("max ADC (ADC)");
1106     fHistDistrMaxIROC->SetYTitle("counts");
1107     fHistDistrMaxIROC->SetTitle(titleMAX);
1108     fHistDistrMaxIROC->Draw("");
1109     fHistDistrMaxOROC->SetLineColor(2);
1110     fHistDistrMaxOROC->Draw("same");
1111     
1112     if(fHistDistrMaxOROC->GetMaximum()> fHistDistrMaxIROC->GetMaximum())  fHistDistrMaxIROC->SetMaximum(fHistDistrMaxOROC->GetMaximum()*1.1);
1113     
1114     TLegend* legio = new TLegend(0.6,0.6,0.8,0.8);
1115     legio->SetFillColor(0);
1116     legio->AddEntry(fHistDistrMaxIROC,"IROC","l");
1117     legio->AddEntry(fHistDistrMaxOROC,"OROC","l");
1118     legio->Draw("same");
1119     
1120       // ADC sum distribution
1121     TCanvas* csum = 0;
1122     csum = (TCanvas*)gROOT->GetListOfCanvases()->FindObject("csum");
1123     if(!csum)  csum = CreateCanvas("csum") ;
1124     csum->cd();
1125     
1126     fHistDistrSumIROC->SetXTitle("sum ADC (ADC)");
1127     fHistDistrSumIROC->SetYTitle("counts");
1128     fHistDistrSumIROC->SetTitle(titleSUM);
1129     fHistDistrSumIROC->Draw("");
1130     fHistDistrSumOROC->SetLineColor(2);
1131     fHistDistrSumOROC->Draw("same");
1132     if(fHistDistrSumOROC->GetMaximum()> fHistDistrSumIROC->GetMaximum())  fHistDistrSumIROC->SetMaximum(fHistDistrSumOROC->GetMaximum()*1.1);
1133     legio->Draw("same");
1134     
1135       // BASELINE MEAN distribution
1136     TCanvas* cbasemean = 0;
1137     cbasemean = (TCanvas*)gROOT->GetListOfCanvases()->FindObject("cbasemean");
1138     if(!cbasemean)  cbasemean = CreateCanvas("cbasemean");
1139     cbasemean->cd();
1140     
1141     fHistDistrBaseMeanIROC = fHistDistrBase2dIROC->ProjectionX("fHistDistrBaseMeanIROC");
1142     fHistDistrBaseMeanIROC->SetXTitle("base mean (ADC)");
1143     fHistDistrBaseMeanIROC->SetYTitle("counts");
1144     fHistDistrBaseMeanIROC->SetTitle(titleMEAN);
1145     fHistDistrBaseMeanIROC->Draw("");
1146     
1147     fHistDistrBaseMeanOROC = fHistDistrBase2dOROC->ProjectionX("fHistDistrBaseMeanOROC");
1148     fHistDistrBaseMeanOROC->SetLineColor(2);
1149     fHistDistrBaseMeanOROC->Draw("same");
1150     if(fHistDistrBaseMeanOROC->GetMaximum()>fHistDistrBaseMeanIROC->GetMaximum())   fHistDistrBaseMeanIROC->SetMaximum(fHistDistrBaseMeanOROC->GetMaximum()*1.1);
1151     legio->Draw("same");
1152     
1153     TCanvas* cbaserms = 0;
1154     cbaserms = (TCanvas*)gROOT->GetListOfCanvases()->FindObject("cbaserms");
1155     if(!cbaserms) cbaserms = CreateCanvas("cbaserms") ;
1156     cbaserms->cd();
1157     
1158       // BASELINE RMS distribution
1159     fHistDistrBaseRmsIROC = fHistDistrBase2dIROC->ProjectionY("fHistDistrBaseRmsIROC");
1160     fHistDistrBaseRmsIROC->SetXTitle("base rms (ADC)");
1161     fHistDistrBaseRmsIROC->SetYTitle("counts");
1162     fHistDistrBaseRmsIROC->SetTitle(titleRMS);
1163     fHistDistrBaseRmsIROC->Draw("");
1164     
1165     fHistDistrBaseRmsOROC = fHistDistrBase2dOROC->ProjectionY("fHistDistrBaseRmsOROC");
1166     fHistDistrBaseRmsOROC->SetLineColor(2);
1167     fHistDistrBaseRmsOROC->Draw("same");
1168     if(fHistDistrBaseRmsOROC->GetMaximum()>fHistDistrBaseRmsIROC->GetMaximum())  fHistDistrBaseRmsIROC->SetMaximum(fHistDistrBaseRmsOROC->GetMaximum()*1.1);
1169     legio->Draw("same");
1170     
1171     cmax->Update();
1172     csum->Update();
1173     cbasemean->Update();
1174     cbaserms->Update();
1175   }
1176   else  if(histos==3)
1177   {
1178       // GLOBAL MAX ADC _________________________________
1179     if(GetProcNextEvent()==1)
1180     {
1181       TCanvas* cglobA =0;
1182       TCanvas* cglobC =0;
1183       
1184       if(!(cglobC=(TCanvas*)gROOT->GetListOfCanvases()->FindObject("SIDE C all"))) cglobC = CreateCanvas("SIDE C all");
1185       if(!(cglobA=(TCanvas*)gROOT->GetListOfCanvases()->FindObject("SIDE A all"))) cglobA = CreateCanvas("SIDE A all");
1186       
1187       Char_t globtitle1[256]; sprintf(globtitle1,"SIDE A Run %05i (EventID %i)",fRunId,fEventNumber);
1188       Char_t globtitle2[256]; sprintf(globtitle2,"SIDE C Run %05i (EventID %i)",fRunId,fEventNumber);
1189       
1190       fHistGlobalMaxA->SetTitle(globtitle1);
1191       fHistGlobalMaxC->SetTitle(globtitle2);
1192       fHistGlobalMaxA->SetXTitle("x/mm");
1193       fHistGlobalMaxA->SetYTitle("y/mm");
1194       fHistGlobalMaxC->SetXTitle("x/mm");
1195       fHistGlobalMaxC->SetYTitle("y/mm");
1196       
1197       if(GetPedestals()==0)     {             fHistGlobalMaxA->SetZTitle("max adc (not baseline sub)");   fHistGlobalMaxC->SetZTitle("max adc (not baseline sub)");  }
1198       else                          {         fHistGlobalMaxA->SetZTitle("max adc ");                     fHistGlobalMaxC->SetZTitle("max adc ");                          }
1199       
1200       fHistGlobalMaxA->SetMinimum(0.01);
1201       fHistGlobalMaxC->SetMinimum(0.01);
1202       
1203       cglobC->cd() ; fHistGlobalMaxC->Draw("COLZ");
1204       cglobA->cd() ; fHistGlobalMaxA->Draw("COLZ");
1205       
1206       Char_t nameom[256];
1207       sprintf(nameom,".x  %s/TPC/AliTPCMonitorExec.C(3)",gSystem->Getenv("ALICE_ROOT"));
1208       
1209       if(fExecGlob==0)
1210       {
1211         if(fVerb)cout << " set exec " << nameom << endl;
1212         cglobC->AddExec("glob",nameom);
1213         cglobA->AddExec("glob",nameom);
1214         fExecGlob = 1;
1215       }
1216       else
1217       {
1218         cglobC->DeleteExec("glob");
1219         cglobA->DeleteExec("glob");
1220         
1221         if(fVerb)  cout << " set exec " << nameom << endl;
1222         cglobC->AddExec("glob",nameom);
1223         cglobA->AddExec("glob",nameom);
1224         
1225       }
1226       cglobC->Update();
1227       cglobA->Update();
1228     }
1229     
1230   }
1231 }
1232
1233
1234
1235 //__________________________________________________________________
1236 void AliTPCMonitor::DrawRMSMap() 
1237 {
1238   // Draw 2Dim rms histos for IROC and OROC
1239   // and set executables for canvases
1240   
1241   TCanvas* crmsoroc =0;
1242   TCanvas* crmsiroc =0;
1243   if(!(crmsoroc = (TCanvas*)gROOT->GetListOfCanvases()->FindObject("crmsoroc")))    crmsoroc    = CreateCanvas("crmsoroc");
1244   if(!(crmsiroc = (TCanvas*)gROOT->GetListOfCanvases()->FindObject("crmsiroc")))    crmsiroc    = CreateCanvas("crmsiroc");
1245   
1246   crmsiroc->cd();  fHistIROCRMS->Draw("COLZ");
1247   crmsoroc->cd();  fHistOROCRMS->Draw("COLZ");
1248   
1249   Char_t carry1[100];  sprintf(carry1,".x %s/TPC/AliTPCMonitorExec.C(1)",gSystem->Getenv("ALICE_ROOT"));
1250   Char_t carry2[100];  sprintf(carry2,".x %s/TPC/AliTPCMonitorExec.C(2)",gSystem->Getenv("ALICE_ROOT"));
1251   
1252   if(fExecPadIrocRms==0)
1253   {
1254     crmsiroc->AddExec("pad",carry1);
1255     crmsiroc->AddExec("row",carry2);
1256     fExecPadIrocRms=1;
1257   }
1258   
1259   if(fExecPadOrocRms==0)
1260   {
1261     crmsoroc->AddExec("pad",carry1);
1262     crmsoroc->AddExec("row",carry2);
1263     fExecPadOrocRms=1;
1264   }
1265   
1266   crmsiroc->Update();
1267   crmsoroc->Update();
1268   
1269   DrawHists(2);
1270   
1271 }
1272
1273 //__________________________________________________________________
1274 void AliTPCMonitor::ExecPad() 
1275 {
1276   
1277   // Executable for Pad
1278   // Show time profile for channel the mouse is pointing at
1279   
1280   Int_t event = gPad->GetEvent();
1281   if (event != 51)   return;
1282   
1283   TObject *select = gPad->GetSelected();
1284   if(!select)    return;
1285   if(!select->InheritsFrom("TH2")) { return;  }
1286   gPad->GetCanvas()->FeedbackMode(kTRUE);
1287   
1288   // get position
1289   Int_t    px        = gPad->GetEventX();
1290   Int_t    py        = gPad->GetEventY();
1291   Float_t  upy       = gPad->AbsPixeltoY(py);
1292   Float_t  upx       = gPad->AbsPixeltoX(px);
1293   Float_t  y         = gPad->PadtoY(upy);
1294   Float_t  x         = gPad->PadtoX(upx);
1295   
1296   Int_t    setrange  = 0;
1297   
1298   TCanvas* cpad      = 0;
1299   //  Char_t   namehist[50];
1300   Char_t   projhist[60];
1301   Char_t   namesel[256];
1302   Char_t   namecanv[256];
1303   
1304   Int_t    xbinmin  = 0;
1305   Int_t    xbinmax  = 0;
1306   Float_t  ybinmin  = 0;
1307   Float_t  ybinmax  = 0;
1308   Int_t    rocid     = 0;
1309   // Check wich Canvas executed the event
1310   TH2S* fHistIndex=0;
1311   sprintf(namesel,select->GetName());
1312   if(strcmp(namesel,"fHistOROC")==0 || strcmp(namesel,"fHistOROCRMS")==0 || strcmp(namesel,"fHistOROCTime")==0 )
1313   {
1314     rocid = 1;
1315     fPadUsedRoc =1;
1316     sprintf(projhist,"ProjectionOROC");
1317     sprintf(namecanv,"coroc_ch");
1318     fHistIndex = fHistOROCIndex;
1319   }
1320   if(strcmp(namesel,"fHistIROC")==0 || strcmp(namesel,"fHistIROCRMS")==0 || strcmp(namesel,"fHistIROCTime")==0 )
1321   {
1322     rocid = 0;
1323     fPadUsedRoc=0;
1324     sprintf(projhist,"ProjectionIROC");
1325     sprintf(namecanv,"ciroc_ch");
1326     fHistIndex = fHistIROCIndex;
1327   }
1328   
1329   // Check if Canvas already existed and get Ranges from former Prjection histogram
1330   if((cpad=(TCanvas*)gROOT->GetListOfCanvases()->FindObject(namecanv)))
1331   {
1332     cpad->cd();
1333     if(gROOT->Get(projhist))
1334     {
1335       setrange = 1;
1336       xbinmin = ((TH1D*)gROOT->Get(projhist))->GetXaxis()->GetFirst();
1337       xbinmax = ((TH1D*)gROOT->Get(projhist))->GetXaxis()->GetLast();
1338       ybinmin = ((TH1D*)gROOT->Get(projhist))->GetMinimum();
1339       ybinmax = ((TH1D*)gROOT->Get(projhist))->GetMaximum();
1340       delete gROOT->Get("legfit");
1341       delete gROOT->Get("fg");
1342     }
1343   }
1344   else
1345   {
1346     cpad = CreateCanvas(namecanv); cpad->cd();
1347   }
1348   
1349   // Get Bin
1350   Int_t testy = fHistIndex->GetYaxis()->FindBin(y);
1351   Int_t testx = fHistIndex->GetXaxis()->FindBin(x);
1352   Int_t binchannel = (Int_t)fHistIndex->GetCellContent(testx,testy);
1353   if(binchannel>30000 || binchannel<0) return;
1354   
1355   TH1D *hp=(TH1D*)gROOT->Get(projhist);
1356   if(!hp)  hp=new TH1D(projhist,projhist,GetTimeBins(),0,GetTimeBins());//delete gROOT->Get(projhist);
1357   // Get Projection
1358 //   TH1D *hp = (TH1D*)(((TH1D*)fHistChannelTime->ProjectionY("hp",binchannel,binchannel))->Clone(projhist));
1359 //   TH1D *hp = fHistChannelTime->ProjectionY(projhist,binchannel,binchannel);
1360 //   hp->Reset();
1361   Int_t nbinsx=fHistChannelTime->GetNbinsX();
1362   for (Int_t itb=0;itb<GetTimeBins();++itb){
1363     hp->GetArray()[itb]=fHistChannelTime->GetArray()[binchannel+itb*(nbinsx+2)];
1364   }
1365 //   printf("hi\n");
1366 //   return;
1367   
1368   // Make title and Pave for channel Info
1369   Char_t title[256];
1370   Int_t npadRow , npad  , nhw , nmax , hwadd;
1371   
1372   hwadd   = (Int_t)fHistChannelTime->GetCellContent(binchannel,0);
1373   fPadUsedHwAddr = hwadd;
1374   
1375   if(rocid==0)npadRow = fMapHand->GetPadRow(hwadd);
1376   else        npadRow = fMapHand->GetPadRow(hwadd)-63;
1377   npad                = fMapHand->GetPad(hwadd);
1378   nhw                 = hwadd;
1379   nmax                = (Int_t)hp->GetMaximum();
1380   
1381   
1382   TPaveText*  legstat = new TPaveText(0.18,0.65,0.3,0.8,"NDC");
1383   
1384   Int_t   connector   =  fMapHand->GetFECconnector(hwadd);
1385   Int_t   fecnr       =  fMapHand->GetFECfromHw(hwadd);
1386   Int_t   fecch       =  fMapHand->GetFECchannel(hwadd);
1387   Int_t   altrochip   =  fMapHand->GetAltro(hwadd);
1388   Int_t   altrochannel=  (fMapHand->GetAltroChannel(hwadd))%16;
1389   Int_t   fecloc      =  fMapHand->U2fGetFECinRCU(fecnr) ;
1390   Int_t   feclocbran  =  fMapHand->U2fGetFECinBranch(fecnr);
1391   Int_t   branch      =  fMapHand->U2fGetBranch(fecnr);
1392   
1393   
1394   Short_t fecget      = (hwadd & fgkHwMaskFEC)   >> 7;
1395   Short_t branchget   = (hwadd & fgkHwMaskBranch)>> 11;
1396   
1397   
1398   Char_t nstat1[100];  Char_t nstat2[100];  Char_t nstat3[100];  Char_t nstat4[100];
1399   Char_t nstat5[100];  Char_t nstat6[100];  Char_t nstat7[100];  Char_t nstat8[100];
1400   
1401   sprintf(nstat1,"Branch (map) \t %i (%i) \n",branchget,branch);
1402   sprintf(nstat2,"Fec in patch \t %i \n",fecloc);
1403   sprintf(nstat8,"Fec in branch (map)\t %i (%i)\n",fecget,feclocbran);
1404   sprintf(nstat7,"Connector  \t %i \n",connector);
1405   sprintf(nstat3,"Fec No.   \t %i \n",fecnr);
1406   sprintf(nstat4,"Fec chan  \t %i \n",fecch);
1407   sprintf(nstat5,"Altro chip\t %i \n",altrochip);
1408   sprintf(nstat6,"Altro chan\t %i \n",altrochannel);
1409   
1410   legstat->AddText(nstat1); legstat->AddText(nstat2);  legstat->AddText(nstat8);  legstat->AddText(nstat7);
1411   legstat->AddText(nstat3); legstat->AddText(nstat4);  legstat->AddText(nstat5);  legstat->AddText(nstat6);
1412   
1413   sprintf(title,"Row=%d Pad=%d Hw =%d maxADC =%d count =%d",npadRow,npad,nhw,nmax,binchannel);
1414   
1415 //   hp->SetName(projhist);
1416   hp->SetTitleSize(0.04);
1417   hp->SetTitle(title);
1418   hp->SetYTitle("ADC");
1419   hp->SetXTitle("Timebin");
1420   hp->GetXaxis()->SetTitleColor(1);
1421   
1422   if(setrange)
1423   {
1424     hp->GetXaxis()->SetRange(xbinmin,xbinmax);
1425     hp->SetMinimum(ybinmin);
1426     hp->SetMaximum(ybinmax);
1427   }
1428   else
1429   {
1430     hp->SetMinimum(0.0);
1431     hp->SetMaximum(1000.0);
1432   }
1433   
1434   cpad->cd();
1435   hp->Draw();
1436   
1437   // Make Fit to peak
1438   if(GetPedestals() && fDisableFit==0)
1439   {
1440     Int_t maxx  =    (Int_t)fHistAddrMaxAdcX->GetBinContent(hwadd);
1441     Float_t max  =  (Float_t)fHistAddrMaxAdc->GetBinContent(hwadd);
1442     Float_t base =  (Float_t)fHistAddrBaseMean->GetBinContent(hwadd);
1443     if(base!=0)
1444     {
1445       if( ((max+base)/base)>1.2)
1446       {
1447         TF1* fg = new TF1("fg",AliTPCMonitor::Gamma4,maxx-5,maxx+5,4);
1448         fg->SetParName(0,"Normalisation");
1449         fg->SetParName(1,"Minimum");
1450         fg->SetParName(2,"Width");
1451         fg->SetParName(3,"Base");
1452         fg->SetParameter(0,max);
1453         fg->SetParameter(1,maxx-2);
1454         fg->SetParameter(2,1.5);
1455         fg->FixParameter(3,0);
1456         fg->SetLineColor(4);
1457         fg->SetLineWidth(1);
1458         hp->Fit("fg","RQ");
1459         
1460         TLegend* legfit = new TLegend(0.6,0.7,0.7,0.8);
1461         legfit->AddEntry("fg","#Gamma 4 fit","l");
1462         legfit->SetFillColor(0);
1463         legfit->SetName("legfit");
1464         legfit->Draw("same");
1465       }
1466     }
1467   }
1468   legstat->SetFillColor(0);
1469   legstat->Draw("same");
1470   cpad->Update();
1471   return;
1472 }
1473
1474 //__________________________________________________________________
1475 void AliTPCMonitor::ExecRow() 
1476 {
1477   
1478   // Executable for Pad
1479   // Show profile of max adc over given pad row
1480   // and 2dim histo  adc(pad-in-row,time bin)
1481   
1482   Int_t event = gPad->GetEvent();
1483   if (event != 61)    return;
1484   gPad->cd();
1485   TObject *select = gPad->GetSelected();
1486   if(!select)    return;
1487   if(!select->InheritsFrom("TH2")) {   return;  }
1488   
1489   Int_t rocid = 0;
1490   //  Char_t namehist[50];
1491   Char_t rowhist[60];
1492   Char_t rowhistsum[60];
1493   Char_t rowhistmax[60];
1494   Char_t rowhistxmax[60];
1495   
1496   sprintf(rowhist,         "hrowtime");
1497   sprintf(rowhistxmax    ,"hxmax");
1498   sprintf(rowhistmax    , "hrowmax");
1499   
1500   // get position
1501   Int_t   px  = gPad->GetEventX();
1502   Int_t   py  = gPad->GetEventY();
1503   Float_t upy = gPad->AbsPixeltoY(py);
1504   Float_t upx = gPad->AbsPixeltoX(px);
1505   Float_t y   = gPad->PadtoY(upy);
1506   Float_t x   = gPad->PadtoX(upx);
1507   
1508   TCanvas*crowtime   = 0;
1509   TCanvas*crowmax    = 0;
1510   TCanvas*cxmax      = 0;
1511   
1512   TH2S*  fHistIndex  = 0;
1513   
1514   // ranges from already existing histos
1515   Int_t    rowtimexmin = 0;
1516   Int_t    rowtimexmax = 0;
1517   Int_t    rowtimeymin = 0;
1518   Int_t    rowtimeymax = 0;
1519   Float_t  rowtimezmin = 0;
1520   Float_t  rowtimezmax = 0;
1521   
1522   Int_t    profrowxmin = 0;
1523   Int_t    profrowxmax = 0;
1524   Double_t profrowymin = 0;
1525   Double_t profrowymax = 0;
1526   
1527   Int_t    profxxmin   = 0;
1528   Int_t    profxxmax   = 0;
1529   Double_t profxymin   = 0;
1530   Double_t profxymax   = 0;
1531   
1532   
1533   Int_t    setrange      = 0;
1534   
1535   
1536   if(     strcmp(select->GetName(),"fHistIROC")==0 || strcmp(select->GetName(),"fHistIROCRMS")==0 ) { fHistIndex = fHistIROCIndex;     rocid =1;   }
1537   else if(strcmp(select->GetName(),"fHistOROC")==0 || strcmp(select->GetName(),"fHistOROCRMS")==0 ) { fHistIndex = fHistOROCIndex;     rocid =2;   }
1538   else                                                                                  { cout << " not implemented for this histo " << endl; return; }
1539   
1540   gPad->GetCanvas()->FeedbackMode(kTRUE);
1541   
1542   
1543   
1544   // check if canvases exist //
1545   crowtime  = (TCanvas*)gROOT->GetListOfCanvases()->FindObject("crowtime");
1546   crowmax   = (TCanvas*)gROOT->GetListOfCanvases()->FindObject("crowmax");
1547   cxmax     = (TCanvas*)gROOT->GetListOfCanvases()->FindObject("cxmax");
1548   
1549   if(!crowtime)   crowtime   = CreateCanvas("crowtime") ;
1550   if(!crowmax)    crowmax    = CreateCanvas("crowmax")  ;
1551   if(!cxmax  )    cxmax      = CreateCanvas("cxmax")    ;
1552   
1553   // check ranges of already existing histos
1554   if(gROOT->Get(rowhist))
1555   {
1556     rowtimexmin  = ((TH2F*)gROOT->Get(rowhist))->GetXaxis()->GetFirst();
1557     rowtimexmax  = ((TH2F*)gROOT->Get(rowhist))->GetXaxis()->GetLast();
1558     rowtimeymin  = ((TH2F*)gROOT->Get(rowhist))->GetYaxis()->GetFirst();
1559     rowtimeymax  = ((TH2F*)gROOT->Get(rowhist))->GetYaxis()->GetLast();
1560     rowtimezmin  = ((TH2F*)gROOT->Get(rowhist))->GetMinimum();
1561     rowtimezmax  = ((TH2F*)gROOT->Get(rowhist))->GetMaximum();
1562     
1563     profrowxmin  = ((TH1F*)gROOT->Get(rowhistmax))->GetXaxis()->GetFirst();
1564     profrowxmax  = ((TH1F*)gROOT->Get(rowhistmax))->GetXaxis()->GetLast();
1565     profrowymin  = ((TH1F*)gROOT->Get(rowhistmax))->GetMinimum();
1566     profrowymax  = ((TH1F*)gROOT->Get(rowhistmax))->GetMaximum();
1567     
1568     profxxmin    = ((TH1F*)gROOT->Get(rowhistxmax))->GetXaxis()->GetFirst();
1569     profxxmax    = ((TH1F*)gROOT->Get(rowhistxmax))->GetXaxis()->GetLast();
1570     profxymin    = ((TH1F*)gROOT->Get(rowhistxmax))->GetMinimum();
1571     profxymax    = ((TH1F*)gROOT->Get(rowhistxmax))->GetMaximum();
1572     
1573     setrange =1;
1574     
1575     delete gROOT->Get(rowhist);
1576     delete gROOT->Get(rowhistmax);
1577     delete gROOT->Get(rowhistsum);
1578     delete gROOT->Get("hxmax");
1579     delete gROOT->Get("legrow");
1580   }
1581   
1582   // get channel for xy bin -> getRow -> getNrows -> getHw for each Pad in Row -> get channel for each hw -> make proj
1583   Int_t testy      = fHistIndex->GetYaxis()->FindBin(y);
1584   Int_t testx      = fHistIndex->GetXaxis()->FindBin(x);
1585   Int_t binchannel = (Int_t)fHistIndex->GetCellContent(testx,testy);
1586   
1587   if(binchannel>30000)    return;
1588   if(binchannel<=0 ) { crowtime->Update() ;    crowmax->Update() ;    return ;  }
1589   
1590   // get hwaddress
1591   Int_t hwadd     = (Int_t)fHistChannelTime->GetCellContent(binchannel,0);
1592   Int_t row       = fMapHand->GetPadRow(hwadd);
1593   Int_t pad       = fMapHand->GetPad(hwadd)   ;
1594   Int_t numofpads =  fMapHand->GetNumofPads(row);
1595   
1596   // create histos
1597   TH2F *hrowtime = new TH2F(rowhist     , ""  ,numofpads,0,numofpads,GetTimeBins(),0.0,GetTimeBins());
1598   TH1F *hrowmax  = new TH1F(rowhistmax , ""  ,numofpads,0,numofpads);
1599   TH1F *hxmax    = new TH1F(rowhistxmax, ""  ,159,0,159      );
1600   
1601   // Row profile ///////////
1602   if(fVerb) cout << " Number of pads " << numofpads << endl;
1603   for(Int_t padnr = 0; padnr<numofpads;padnr++)
1604   {
1605     Int_t addrinrow = fMapHand->GetPadAddInRow(row,padnr );
1606     Int_t channel   = (Int_t)fHistAddrMapIndex->GetBinContent(addrinrow);
1607     if(channel==-1) continue;
1608     
1609     hrowmax->SetBinContent(padnr+1,fHistAddrMaxAdc->GetBinContent(addrinrow));
1610     TH1D *hp = fHistChannelTime->ProjectionY("hp",channel,channel);
1611     for(Int_t time = 0;time<GetTimeBins();time++) {
1612       
1613       Float_t val = hp->GetBinContent(time);
1614       hrowtime->SetCellContent(padnr+1,time+1,val);
1615     }
1616   }
1617   
1618   // X profile  /////////////
1619   Double_t xval  = 0.0;
1620   Double_t yval  = 0.0;
1621   GetXY(xval,yval,numofpads,row,pad);
1622   
1623   Int_t padnr = 0;
1624   Int_t hw    = 0;
1625   for(Int_t nrow = 0; nrow<159; nrow++)
1626   {
1627     padnr = GetPadAtX(xval,nrow);
1628     if(padnr>=0)
1629     {
1630       hw = fMapHand->GetPadAddInRow(nrow,padnr);
1631       if(fPadMapHw[hw]==-1){ continue                                                      ; }
1632       else                { hxmax->SetBinContent(nrow+1,fHistAddrMaxAdc->GetBinContent(hw))   ; }
1633     }
1634   }
1635   
1636   cxmax->cd();
1637   Char_t hxtitle[50] ; sprintf(hxtitle,"max adc in pads at x=%5.1f mm",xval);
1638   hxmax->SetTitle(hxtitle);
1639   hxmax->SetXTitle("row");
1640   if(!GetPedestals()) hxmax->SetYTitle("max adc (baseline sub.)");
1641   else                hxmax->SetYTitle("max adc ");
1642   hxmax->SetMinimum(0.01);
1643   hxmax->Draw("l");
1644   
1645   if(setrange)
1646   {
1647     hxmax->GetXaxis()->SetRange(profxxmin,profxxmax);
1648     hxmax->SetMinimum(profxymin);
1649     hxmax->SetMaximum(profxymax);
1650   }
1651   
1652   cxmax->Update();
1653   
1654   crowtime->cd();
1655   Char_t title[256];
1656   Char_t titlemax[256];
1657   if(rocid==1) {sprintf(title,"%s Row=%d",((TH2*)select)->GetTitle(),row)   ;    sprintf(titlemax,"IROC  max/sum Row=%d",row   );}
1658   else         {sprintf(title,"%s Row=%d",((TH2*)select)->GetTitle(),row-63);    sprintf(titlemax,"OROC  max/sum Row=%d",row-63);}
1659   if(fVerb) cout << " set name " << endl;
1660   
1661   
1662   // row vs time
1663   crowtime->cd();
1664   hrowtime->SetTitleSize(0.04);
1665   hrowtime->SetTitle(title);
1666   hrowtime->SetYTitle("timbin");
1667   hrowtime->SetXTitle("pad in row");
1668   hrowtime->SetZTitle("signal (ADC)");
1669   
1670   hrowtime->GetXaxis()->SetTitleColor(1);
1671   hrowtime->SetMaximum(1000.0);
1672   hrowtime->SetMinimum(0.0);
1673   
1674   if(setrange)
1675   {
1676     hrowtime->GetXaxis()->SetRange(rowtimexmin,rowtimexmax);
1677     hrowtime->GetYaxis()->SetRange(rowtimeymin,rowtimeymax);
1678     hrowtime->SetMinimum(rowtimezmin);
1679     hrowtime->SetMaximum(rowtimezmax);
1680   }
1681   
1682   hrowtime->Draw("COLZ");
1683   crowtime->UseCurrentStyle();
1684   crowtime->Update();
1685   
1686   // max and sum /////////////////////////
1687   crowmax->cd();
1688   if(setrange) {
1689     hrowmax->GetXaxis()->SetRange(profrowxmin,profrowxmax);
1690     hrowmax->SetMinimum(profrowymin);
1691     hrowmax->SetMaximum(profrowymax);
1692   }
1693   hrowmax->SetTitleSize(0.04);
1694   hrowmax->SetTitle(title);
1695   hrowmax->SetYTitle("max adc");
1696   hrowmax->SetXTitle("pad in row");
1697   hrowmax->GetXaxis()->SetTitleColor(1);
1698   
1699   hrowmax->SetLineColor(2);
1700   hrowmax->Draw("l");
1701   crowmax->Update();
1702   
1703   return;
1704 }
1705
1706 //__________________________________________________________________
1707 void AliTPCMonitor::Write10bitChannel()
1708 {
1709   
1710   // Write 10 bit words form histogram for active(last pointed)  channel
1711   
1712   if(fPadUsedHwAddr==-1){ AliWarning(" No active pad "); return ;}
1713   
1714   Int_t  pad     = (Int_t)fMapHand->GetPad(   fPadUsedHwAddr);
1715   Int_t  row     = (Int_t)fMapHand->GetPadRow(fPadUsedHwAddr);
1716   Int_t  channel = (Int_t)fPadMapHw[fPadUsedHwAddr];
1717   
1718   Char_t filenameroot[256];
1719   Char_t filenamedat[256];
1720   Char_t projhist[256];
1721   
1722   if(fPadUsedRoc==1) { sprintf(projhist,"ProjectionOROC"); }
1723   if(fPadUsedRoc==0) { sprintf(projhist,"ProjectionIROC"); }
1724   
1725   sprintf(filenamedat, "Channel_Run%05i_EventID_%i_Pad_%i_Row_%i.dat"      ,fRunId,fEventNumber,pad,row);
1726   sprintf(filenameroot,"Channel_Run%05i_EventID_%i_Pad_%i_Row_%i.root"     ,fRunId,fEventNumber,pad,row);
1727   
1728   TH1D* hpr = 0;
1729   if((hpr=(TH1D*)gROOT->Get(projhist)))
1730   {
1731       // root histo
1732     TFile f(filenameroot,"recreate");
1733     hpr->Write();
1734     f.Close();
1735     
1736       // raw singal
1737     ofstream datout(filenamedat,ios::out);
1738     datout <<"Timebin \t ADC value " << endl;
1739     for(Int_t i = 1; i <GetTimeBins(); i++)
1740     {
1741       datout << i << " \t \t " << fPad[channel][i] << endl;
1742     }
1743     datout.close();
1744   }
1745   else
1746   {
1747     AliWarning("No projection histo found ");
1748   }
1749 }
1750
1751 //__________________________________________________________________
1752 void AliTPCMonitor::ExecTransform() 
1753 {
1754   
1755   // Make Fast Fourier Transformation for active pad
1756   // fft is only performed for a data sample of size 2^n
1757   // reduce window according to largest  power of 2 which is smaller than the viewing  range
1758   
1759   Char_t namecanv[256];
1760   Char_t namecanv2[256];
1761   Char_t projhist[256];
1762   Char_t namehtrimag[256];
1763   Char_t namehtrreal[256];
1764   Char_t namehtrmag[256];
1765   
1766   if(fPadUsedRoc==1) {    sprintf(namecanv,"coroc_ch_trans") ; sprintf(namecanv2,"coroc_ch_trans2") ;   sprintf(projhist,"ProjectionOROC");  }
1767   if(fPadUsedRoc==0) {    sprintf(namecanv,"ciroc_ch_trans") ; sprintf(namecanv2,"ciroc_ch_trans2") ;  sprintf(projhist,"ProjectionIROC");  }
1768   
1769   TH1D*  hproj = 0;
1770   
1771   if((TH1D*)gROOT->Get(projhist)==0){AliWarning("Proj histo does not exist \n Move mouse over 2d histo choose channel \n and drag mouse form histo again!");  return ;}
1772   else      hproj = (TH1D*)gROOT->Get(projhist) ;
1773   
1774   
1775   if(fPadUsedRoc==1) {  sprintf(namehtrimag,"htransimagfreq_oroc");    sprintf(namehtrreal,"htransrealfreq_oroc"); sprintf(namehtrmag,"htransmagfreq_oroc");  }
1776   else               {  sprintf(namehtrimag,"htransimagfreq_iroc");    sprintf(namehtrreal,"htransrealfreq_iroc"); sprintf(namehtrmag,"htransmagfreq_iroc");  }
1777   
1778   if( gROOT->Get(namehtrimag))  delete  gROOT->Get(namehtrimag);
1779   if( gROOT->Get(namehtrreal))  delete  gROOT->Get(namehtrreal);
1780   if( gROOT->Get(namehtrmag))  delete  gROOT->Get(namehtrmag);
1781   
1782   TCanvas *ctrans = 0;
1783   if(!(ctrans = (TCanvas*)gROOT->GetListOfCanvases()->FindObject(namecanv)))
1784   {
1785     ctrans = CreateCanvas(namecanv);
1786     ctrans->Divide(1,2);
1787   }
1788   TCanvas *ctrans2 = 0;
1789   if(!(ctrans2 = (TCanvas*)gROOT->GetListOfCanvases()->FindObject(namecanv2)))
1790   {
1791     ctrans2 = CreateCanvas(namecanv2);
1792 //      ctrans2->Divide(1,2);
1793   }
1794   
1795   Int_t binfirst  =  hproj->GetXaxis()->GetFirst();
1796   Int_t binlast   =  hproj->GetXaxis()->GetLast();
1797   Int_t bins       =  binlast -binfirst +1;
1798   
1799   Int_t power = 0;
1800   for(Int_t pot = 0; pot<=10 ; pot++)
1801   {
1802     Int_t comp =  (Int_t)TMath::Power(2,pot);
1803     if(bins>=comp)power = pot;
1804   }
1805   
1806   bins = (Int_t)TMath::Power(2,power);
1807   
1808   // sampling frequency ;
1809   Double_t  deltat = 1.0/(Float_t)GetSamplingFrequency();
1810   
1811   // output histo
1812   TH1D* htransrealfreq = new TH1D(namehtrreal,namehtrreal,10000,-1/(2*deltat),1/(2*deltat));
1813   TH1D* htransimagfreq = new TH1D(namehtrimag,namehtrimag,10000,-1/(2*deltat),1/(2*deltat));
1814   TH1D* htransmag      = new TH1D(namehtrmag,namehtrmag,10000,-1/(2*deltat),1/(2*deltat));
1815   
1816   Char_t titlereal[256];
1817   Char_t titleimag[256];
1818   Char_t titlemag[256];
1819   if(fPadUsedRoc==1) {    sprintf(titlereal,"OROC DFT real part");  sprintf(titleimag,"OROC DFT imag part");  sprintf(titlemag,"OROC DFT magnitude");  }
1820   else {                  sprintf(titlereal,"IROC DFT real part");  sprintf(titleimag,"IROC DFT imag part");  sprintf(titlemag,"IROC DFT magnitude");  }
1821   
1822   htransrealfreq->SetTitle(titlereal);  htransrealfreq->SetXTitle("f/hz");  htransrealfreq->SetYTitle("z_{real}(f)");
1823   htransimagfreq->SetTitle(titleimag);  htransimagfreq->SetXTitle("f/hz");  htransimagfreq->SetYTitle("z_{imag}(f)");
1824   htransmag->SetTitle(titlemag);  htransmag->SetXTitle("f/hz");  htransmag->SetYTitle("mag(f)");
1825   
1826   // create complex packed data array
1827   const Int_t kdatasiz = 2*bins;
1828   Double_t*  data = new Double_t[kdatasiz];
1829   for(Int_t i=0;i<2*bins;i++)  { data[i]   =  0.0;}
1830   for(Int_t i=0;i<bins;i++)    { data[2*i] = (Double_t)hproj->GetBinContent(binfirst+i); }
1831   
1832   // make fourier transformation
1833   AliTPCMonitorFFT* four = new AliTPCMonitorFFT();
1834   four->ComplexRadix2ForwardWrap(data,1,bins);
1835   
1836   // write output  and fill histos forward
1837   Double_t freq =  0.0;
1838   for(Int_t i=0;i<2*bins;i++)
1839   {
1840     if(i<bins)
1841     {
1842       if(i<(bins/2))  { freq = i/(bins*deltat)            ; }
1843       else            { freq = -1*((bins-i)/(bins*deltat)); }
1844       htransrealfreq->Fill( freq,data[2*i]  );
1845       htransimagfreq->Fill( freq,data[2*i+1]);
1846       htransmag->Fill( freq, TMath::Sqrt(data[2*i]*data[2*i]+data[2*i+1]*data[2*i+1]) );
1847     }
1848   }
1849   
1850   ctrans->cd(1);
1851   htransrealfreq->Draw();
1852   ctrans->cd(2);
1853   htransimagfreq->Draw();
1854   ctrans->Update();
1855   ctrans2->cd();
1856   htransmag->Draw();
1857   ctrans2->Update();
1858   delete four;
1859   delete data;
1860 }
1861
1862 //__________________________________________________________________
1863 void AliTPCMonitor::ShowSel(Int_t* compval)               
1864 {
1865   
1866   // Show only selected components
1867   // First restore original histogram from clone
1868   // Than remove all not matching pads form histos
1869   
1870   Int_t   connector   =  0;
1871   Int_t   fecnr       =  0;
1872   Int_t   altrochip   =  0;
1873   Int_t   feclocbran  =  0;
1874   Int_t   branch      =  0;
1875   Short_t rcuget      =  0;
1876   Int_t   emptyI      =  1;
1877   Int_t   index       = -1;
1878   Int_t   hwadd       =  0;
1879   
1880   Float_t maxiroc     =  fHistIROCClone->GetMaximum();
1881   Float_t maxoroc     =  fHistOROCClone->GetMaximum();
1882   
1883   
1884   //  restore original histos
1885   for(Int_t row = 0; row<fkNRowsIroc;  row++)
1886   {
1887     for(Int_t pad = 0; pad<fkNPadsIroc;  pad++)
1888     {
1889       index = (Int_t)fHistIROCIndex->GetCellContent(row+1,pad+1);
1890       if(index==-1)continue;
1891       else  fHistIROC->SetCellContent(row+1,pad+1,fHistIROCClone->GetCellContent(row+1,pad+1));
1892     }
1893   }
1894   for(Int_t row = 0; row<fkNRowsOroc;  row++)
1895   {
1896     for(Int_t pad = 0; pad<fkNPadsOroc;  pad++)
1897     {
1898       index = (Int_t)fHistOROCIndex->GetCellContent(row+1,pad+1);
1899       if(index==-1)continue;
1900       else    fHistOROC->SetCellContent(row+1,pad+1,fHistOROCClone->GetCellContent(row+1,pad+1));
1901     }
1902   }
1903   
1904   
1905   // remove not matching entries from fHistIROC/fHistOROC
1906   
1907   TH2F* fHist       =0;
1908   TH2S* fHistIndex  =0;
1909   Int_t npads       =0;
1910   Int_t subrows     =0;
1911   
1912   for(Int_t row = 0; row< (fkNRowsIroc + fkNRowsOroc);  row++)
1913   {
1914     if(row<fkNRowsIroc) {  fHist=fHistIROC ; fHistIndex = fHistIROCIndex; npads = fkNPadsIroc; subrows =0         ;}
1915     else                {  fHist=fHistOROC ; fHistIndex = fHistOROCIndex; npads = fkNPadsOroc; subrows =fkNRowsIroc;}
1916     
1917     for(Int_t pad = 0; pad<npads;  pad++)
1918     {
1919       index    = (Int_t)fHistIndex->GetCellContent(row -subrows +1,pad+1);
1920       if(index==-1)  continue ;
1921       hwadd    = (Int_t)fHistChannelTime->GetCellContent(index,0);
1922       
1923     // global fecnr
1924       fecnr     =  fMapHand->GetFECfromHw(hwadd);
1925       if(compval[0]!=-1 && fecnr!=compval[0])      {      fHist->SetCellContent(row-subrows+1,pad+1,0);           continue;     }
1926       
1927     // rcu
1928       rcuget      = (hwadd & fgkHwMaskRCU)>> 12;
1929       if(compval[1]!=-1 && rcuget!=compval[1])     {      fHist->SetCellContent(row-subrows+1,pad+1,0);           continue;     }
1930       
1931     // branch
1932       branch    =  fMapHand->U2fGetBranch(fecnr) ;
1933       if(compval[2]!=-1 && branch!=compval[2])     {      fHist->SetCellContent(row-subrows+1,pad+1,0);           continue;     }
1934       
1935     // local fec
1936       feclocbran=   fMapHand->U2fGetFECinBranch(fecnr) ;
1937       if(compval[3]!=-1 && feclocbran!=compval[3]) {      fHist->SetCellContent(row-subrows+1,pad+1,0);           continue;     }
1938       
1939     // connector
1940       connector =  fMapHand->GetFECconnector(hwadd);
1941       if(compval[4]!=-1 && connector!=compval[4])  {      fHist->SetCellContent(row-subrows+1,pad+1,0);           continue;     }
1942       
1943     // Altro chip
1944       altrochip     =  fMapHand->GetAltro(hwadd);
1945       if(compval[5]!=-1 && altrochip!=compval[5])      {          fHist->SetCellContent(row-subrows+1,pad+1,0);           continue;     }
1946       emptyI =0;
1947     }
1948   }
1949   
1950   TCanvas* c1 =0;
1951   TCanvas* c2 =0;
1952   if(gROOT->GetListOfCanvases()->FindObject("ciroc"))
1953   {
1954     c1 = (TCanvas*)gROOT->GetListOfCanvases()->FindObject("ciroc");
1955     c1->cd() ;
1956     fHistIROC->Draw("COLZ");
1957     fHistIROC->SetMaximum(maxiroc);
1958     fHistIROC->SetMinimum(0.0);
1959     c1->Update();
1960   }
1961   if(gROOT->GetListOfCanvases()->FindObject("coroc"))
1962   {
1963     c2 =  (TCanvas*)gROOT->GetListOfCanvases()->FindObject("coroc");
1964     c2->cd() ;
1965     fHistOROC->Draw("COLZ");
1966     fHistOROC->SetMaximum(maxoroc);
1967     fHistOROC->SetMinimum(0.0);
1968     c2->Update();
1969   }
1970   return ;
1971 }
1972
1973 //__________________________________________________________________
1974 void AliTPCMonitor::ResizeCanv() 
1975 {
1976   // Resize canvases and delete some of them
1977   
1978   Char_t carry1[100];
1979   sprintf(carry1,".x %s/TPC/AliTPCMonitorExec.C(1)",gSystem->Getenv("ALICE_ROOT"));
1980   Char_t carry3[100];
1981   sprintf(carry3,".x %s/TPC/AliTPCMonitorExec.C(2)",gSystem->Getenv("ALICE_ROOT"));
1982   if(fVerb) cout <<  " canv 1 " << endl;
1983   
1984   if(gROOT->GetListOfCanvases()->FindObject(        "coroc_ch")) {  delete gROOT->GetListOfCanvases()->FindObject("coroc_ch") ; }
1985   if(gROOT->GetListOfCanvases()->FindObject(        "ciroc_ch")) {  delete gROOT->GetListOfCanvases()->FindObject("ciroc_ch") ; }
1986   
1987   // for 2dim plots delete create and draw again
1988   if(gROOT->GetListOfCanvases()->FindObject("ciroc"))
1989   {
1990     delete  gROOT->GetListOfCanvases()->FindObject("ciroc");
1991     TCanvas* ciroc = CreateCanvas("ciroc");
1992     ciroc->cd();
1993     fHistIROC->Draw("COLZ");
1994     ciroc->AddExec("pad",carry1);
1995     ciroc->AddExec("row",carry3);
1996     fExecPlaneMax=1;
1997     ciroc->Update();
1998   }
1999   // for 2dim plots delete create and draw again
2000   if(gROOT->GetListOfCanvases()->FindObject("coroc"))
2001   {
2002     delete gROOT->GetListOfCanvases()->FindObject("coroc");
2003     TCanvas* coroc = CreateCanvas("coroc");
2004     coroc->cd();
2005     fHistOROC->Draw("COLZ");
2006     
2007     coroc->AddExec("pad",carry1);
2008     coroc->AddExec("row",carry3);
2009     coroc->Update();
2010     fExecPlaneMax=1;
2011   }
2012   
2013   if(gROOT->GetListOfCanvases()->FindObject(       "cbasemean")) {    delete gROOT->GetListOfCanvases()->FindObject("cbasemean"); }
2014   if(gROOT->GetListOfCanvases()->FindObject(           "cbase")) {    delete gROOT->GetListOfCanvases()->FindObject("cbase");}
2015   if(gROOT->GetListOfCanvases()->FindObject(        "cbaserms")) {    delete gROOT->GetListOfCanvases()->FindObject("cbaserms");  }
2016   if(gROOT->GetListOfCanvases()->FindObject(            "cmax")) {    delete gROOT->GetListOfCanvases()->FindObject("cmax");      }
2017   if(gROOT->GetListOfCanvases()->FindObject(            "csum")) {    delete gROOT->GetListOfCanvases()->FindObject("csum");      }
2018   if(gROOT->GetListOfCanvases()->FindObject(  "ciroc_ch_trans")) {    delete gROOT->GetListOfCanvases()->FindObject("ciroc_ch_trans");}
2019   if(gROOT->GetListOfCanvases()->FindObject(  "coroc_ch_trans")) {    delete gROOT->GetListOfCanvases()->FindObject("coroc_ch_trans");}
2020   if(gROOT->GetListOfCanvases()->FindObject(       "crowtime")) {    delete gROOT->GetListOfCanvases()->FindObject("crowtime"); }
2021   if(gROOT->GetListOfCanvases()->FindObject(        "crowmax")) {    delete gROOT->GetListOfCanvases()->FindObject("crowmax");  }
2022   if(gROOT->GetListOfCanvases()->FindObject(          "cxmax")) {    delete gROOT->GetListOfCanvases()->FindObject("cxmax");    }
2023   if(gROOT->GetListOfCanvases()->FindObject(        "crmsoroc")) {    delete gROOT->GetListOfCanvases()->FindObject("crmsoroc");         fExecPadOrocRms = 0;   }
2024   if(gROOT->GetListOfCanvases()->FindObject(        "crmsiroc")) {    delete gROOT->GetListOfCanvases()->FindObject("crmsiroc");         fExecPadIrocRms = 0;   }
2025   if(gROOT->GetListOfCanvases()->FindObject(        "crowmax")) {    delete gROOT->GetListOfCanvases()->FindObject("crowmax");  }
2026   if(gROOT->GetListOfCanvases()->FindObject(        "crowmax")) {    delete gROOT->GetListOfCanvases()->FindObject("crowmax");  }
2027   
2028 }
2029
2030
2031
2032
2033 //__________________________________________________________________
2034 Int_t AliTPCMonitor::ExecProcess() 
2035 {
2036   // Executable for global Histogram
2037   // Will be called from /TPC/AliTPCMonitorExec.C(3)
2038   // Call ProcessEvent for same event and sector pointed at
2039   
2040   Int_t side   = 0;
2041   Int_t sector = 0;
2042   
2043   Int_t event = gPad->GetEvent();
2044   if(event != 61)  return -1;
2045   
2046   TObject *select = gPad->GetSelected();
2047   if(!select)  return -1;
2048   if(!select->InheritsFrom("TH2")) {gPad->SetUniqueID(0);    return -1;  }
2049   if(       strcmp(select->GetName(),"hglobal" )==0 || ( strcmp(select->GetName(),"SIDE A" )==0) ) side = 0;
2050   else  if( strcmp(select->GetName(),"hglobal2")==0 || ( strcmp(select->GetName(),"SIDE C" )==0) ) side = 1;
2051   
2052   // get position
2053   Int_t   px    = gPad->GetEventX();
2054   Int_t   py    = gPad->GetEventY();
2055   Float_t upy   = gPad->AbsPixeltoY(py);
2056   Float_t upx   = gPad->AbsPixeltoX(px);
2057   Float_t y     = gPad->PadtoY(upy);
2058   Float_t x     = gPad->PadtoX(upx);
2059   
2060   Int_t testy = ((TH2*)select)->GetYaxis()->FindBin(y);
2061   Int_t testx = ((TH2*)select)->GetXaxis()->FindBin(x);
2062   if(((TH2*)select)->GetCellContent(testx,testy)==0) return -1 ;
2063   
2064   Float_t alpha = 0.0;
2065   if(x>0.0 && y > 0.0)    alpha = TMath::Abs(TMath::ATan(TMath::Abs(x/y)));
2066   if(x>0.0 && y < 0.0)    alpha = TMath::Abs(TMath::ATan(TMath::Abs(y/x)));
2067   if(x<0.0 && y < 0.0)    alpha = TMath::Abs(TMath::ATan(TMath::Abs(x/y)));
2068   if(x<0.0 && y > 0.0)    alpha = TMath::Abs(TMath::ATan(TMath::Abs(y/x)));
2069   
2070   if(x>0.0 && y < 0.0)    alpha += (    TMath::Pi()/2);
2071   if(x<0.0 && y < 0.0)    alpha += (    TMath::Pi());
2072   if(x<0.0 && y > 0.0)    alpha += (1.5*TMath::Pi());
2073   
2074   sector =   (Int_t)(alpha/(2*TMath::Pi()/18.0));
2075   if(alpha> (sector+0.5)*(2*TMath::Pi()/18.0))  sector+=1;
2076   
2077   if(sector==18 && side ==0 ) {
2078     AliWarning("There was a wromg assignment of sector 0 with sector 18. Check sectors");
2079     sector =0;
2080   }
2081   
2082   sector = (18-sector +4)%18;
2083   SetLastSector(sector+ side*18);
2084   SetProcNextEvent(0);
2085   
2086   if(fVerb) cout << "AliTPCMonitor::ExecProcess()  next side          " <<   side    << " next sector        " <<    sector  << endl;
2087   
2088   return (Int_t)ProcessEvent();
2089   
2090 }
2091
2092 //__________________________________________________________________
2093 Int_t AliTPCMonitor::GetRCUPatch(Int_t runid, Int_t eqid) const
2094 {
2095   
2096   // Return RCU patch index for given equipment id eqid
2097   Int_t patch = 0;
2098   //if(runid>=704)
2099   if ( eqid<768 || eqid>983 ) return 0; //no TPC eqid
2100   if(runid>=0)
2101   {
2102     if(eqid>=1000) return 0;
2103     patch = fMapEqidsRcu[eqid] ;
2104   }
2105   else
2106   {
2107     if(eqid==408) {patch =  13*6+4 +0;  }
2108     if(eqid==409) {patch =  13*6+5 +0;  }
2109     if(eqid==509) {patch =  13*6+0 +0;  }
2110     if(eqid==512) {patch =  13*6+3 +0;  }
2111     if(eqid==513) {patch =  13*6+1 +0;  }
2112     if(eqid==517) {patch =  13*6+2 +0;  }
2113     
2114     if(eqid==404) {patch =  4*6+5 +0;   }
2115     if(eqid==504) {patch =  4*6+4 +0;   }
2116     if(eqid==407) {patch =  4*6+3 +0;   }
2117     if(eqid==503) {patch =  4*6+2 +0;   }
2118     if(eqid==508) {patch =  4*6+1 +0;   }
2119     if(eqid==506) {patch =  4*6+0 +0;   }
2120   }
2121   return patch;
2122 }
2123
2124 //__________________________________________________________________
2125 void AliTPCMonitor::DumpHeader(AliRawReader * reader) const
2126 {
2127   // Dump Event header for  format
2128   
2129   cout << "EventHeader     : fReader->GetEquipmentSize()            :" << reader->GetEquipmentSize()        << endl;
2130   cout << "EventHeader     : fReader->GetType()                     :" << reader->GetType()                 << endl;
2131   cout << "EventHeader     : fReader->GetRunNumber()                :" << reader->GetRunNumber()            << endl;
2132   cout << "EventHeader     : fReader->GetEventId()                  :" << *(reader->GetEventId())           << endl;
2133   cout << "EventHeader     : fReader->GetLDCId()                    :" << reader->GetLDCId()                << endl;
2134   cout << "EventHeader     : fReader->GetGDCId()                    :" << reader->GetGDCId()                << endl;
2135 }
2136
2137
2138 //__________________________________________________________________
2139 Double_t AliTPCMonitor::Gamma4(Double_t* x, Double_t* par) {
2140   
2141   // Gamma4 function used to fit signals
2142   // Defined in sections: diverging branch set to 0
2143   
2144   Double_t val  = 0.0;
2145   if(x[0] > par[1])
2146     val = par[0]*exp(4.0)* pow((x[0]-par[1])/par[2],4)*exp(-4.0*(x[0]-par[1])/par[2])+ par[3];
2147   else
2148     val = 0;
2149   return val;
2150 }
2151
2152 //__________________________________________________________________
2153 TCanvas* AliTPCMonitor::CreateCanvas(const Char_t* name)
2154 {
2155   // Create Canvases
2156   
2157   TCanvas* canv =0;
2158   
2159   Int_t xoffset  = GetCanvasXOffset();
2160   Int_t xsize    = GetCanvasXSize();
2161   Int_t ysize    = GetCanvasYSize();
2162   Int_t xspace   = GetCanvasXSpace();
2163   Int_t yspace   = GetCanvasYSpace();
2164   
2165   // ROC 2dim max distribution
2166   if(     strcmp(name,"coroc"         )==0) {    canv   = new TCanvas("coroc"         ,"coroc"      ,                   -1+xoffset,(Int_t)(yspace+0.5*ysize) ,(Int_t)(1.5*xsize),(Int_t)(1.5*ysize)); return canv;  }
2167   else if(strcmp(name,"ciroc"         )==0) {    canv   = new TCanvas("ciroc"         ,"ciroc"      ,                   -1+xoffset,                     0    ,(Int_t)(1.5*xsize),(Int_t)(1.5*ysize)); return canv;  }
2168   // ROC  2dim rms distribution
2169   else if(strcmp(name,"crmsoroc"      )==0) {    canv   = new TCanvas("crmsoroc"      ,"crmsoroc"   ,                   -1+xoffset,(Int_t)(yspace+0.5*ysize) ,(Int_t)(1.5*xsize),(Int_t)(1.5*ysize)); return canv;  }
2170   else if(strcmp(name,"crmsiroc"      )==0) {    canv   = new TCanvas("crmsiroc"      ,"crmsiroc"   ,                   -1+xoffset,                        0 ,(Int_t)(1.5*xsize),(Int_t)(1.5*ysize)); return canv;  }
2171   // Global ADC max Histos
2172   else if(strcmp(name,"SIDE C all"    )==0) {    canv   = new TCanvas("SIDE C all"    ,"SIDE C all" ,   (Int_t)(3*xspace+ xoffset),(Int_t)(yspace+0.5*ysize) ,(Int_t)(1.5*xsize),(Int_t)(1.5*ysize)); return canv;  }
2173   else if(strcmp(name,"SIDE A all"    )==0) {    canv   = new TCanvas("SIDE A all"    ,"SIDE A all" ,   (Int_t)(3*xspace+ xoffset),                        0 ,(Int_t)(1.5*xsize),(Int_t)(1.5*ysize)); return canv;  }
2174   // 1 dim max sum basekine distribution
2175   else if(strcmp(name,"cmax"          )==0) {    canv   = new TCanvas("cmax"          ,"cmax"       ,                   -1+xoffset,                 3*yspace ,(Int_t)(1.0*xsize),(Int_t)(1.0*ysize)); return canv;  }
2176   else if(strcmp(name,"csum"          )==0) {    canv   = new TCanvas("csum"          ,"csum"       ,               xspace+xoffset,                 3*yspace ,(Int_t)(1.0*xsize),(Int_t)(1.0*ysize)); return canv;  }
2177   else if(strcmp(name,"cbasemean"     )==0) {    canv   = new TCanvas("cbasemean"     ,"cbasemean"  ,             2*xspace+xoffset,                 3*yspace ,(Int_t)(1.0*xsize),(Int_t)(1.0*ysize)); return canv;  }
2178   else if(strcmp(name,"cbaserms"      )==0) {    canv   = new TCanvas("cbaserms"      ,"cbaserms"   ,             3*xspace+xoffset,                 3*yspace ,(Int_t)(1.0*xsize),(Int_t)(1.0*ysize)); return canv;  }
2179   // Projections of single channel
2180   else if(strcmp(name,"coroc_ch"      )==0) {    canv   = new TCanvas("coroc_ch"      ,"coroc_ch"   ,   (Int_t)(1.5*xspace+xoffset),(Int_t)(yspace+0.5*ysize),(Int_t)(1.5*xsize),(Int_t)(1.5*ysize)); return canv;  }
2181   else if(strcmp(name,"ciroc_ch"      )==0) {    canv   = new TCanvas("ciroc_ch"      ,"ciroc_ch"   ,   (Int_t)(1.5*xspace+xoffset),                       0 ,(Int_t)(1.5*xsize),(Int_t)(1.5*ysize)); return canv;  }
2182   // FFT for single channel
2183   else if(strcmp(name,"coroc_ch_trans")==0) {    canv   = new TCanvas("coroc_ch_trans","coroc_ch_trans",(Int_t)(3.0*xspace+xoffset),(Int_t)(yspace+0.5*ysize),(Int_t)(1.5*xsize),(Int_t)(1.5*ysize)); return canv;  }
2184   else if(strcmp(name,"ciroc_ch_trans")==0) {    canv   = new TCanvas("ciroc_ch_trans","ciroc_ch_trans",(Int_t)(3.0*xspace+xoffset),                       0 ,(Int_t)(1.5*xsize),(Int_t)(1.5*ysize)); return canv;  }
2185   else if(strcmp(name,"coroc_ch_trans2")==0) {    canv   = new TCanvas("coroc_ch_trans2","coroc_ch_trans2",(Int_t)(3.0*xspace+xoffset),(Int_t)(yspace+0.5*ysize),(Int_t)(1.5*xsize),(Int_t)(1.5*ysize)); return canv;  }
2186   else if(strcmp(name,"ciroc_ch_trans2")==0) {    canv   = new TCanvas("ciroc_ch_trans2","ciroc_ch_trans2",(Int_t)(3.0*xspace+xoffset),                       0 ,(Int_t)(1.5*xsize),(Int_t)(1.5*ysize)); return canv;  }
2187   // row profile histograms
2188   else if(strcmp(name,"crowtime"     )==0) {    canv   = new TCanvas("crowtime"     ,"crowtime"  ,              1*xspace+xoffset,         2*yspace +ysize ,(Int_t)(1.0*xsize),(Int_t)(1.0*ysize)); return canv;  }
2189   else if(strcmp(name,"crowmax"      )==0) {    canv   = new TCanvas("crowmax"      ,"crowmax"   ,              2*xspace+xoffset,         2*yspace +ysize ,(Int_t)(1.0*xsize),(Int_t)(1.0*ysize)); return canv;  }
2190   else if(strcmp(name,"cxmax"        )==0) {    canv   = new TCanvas("cxmax"        ,"cxmax"     ,              3*xspace+xoffset,         2*yspace +ysize ,(Int_t)(1.0*xsize),(Int_t)(1.0*ysize)); return canv;  }
2191   else                                      {    cout   << " Warning Canvas name unknown "  << endl;                                                                                                  return 0   ;  }
2192 }
2193
2194 //__________________________________________________________________
2195 void AliTPCMonitor::WriteHistos()  
2196 {
2197   // Writes all available histograms to a file in current directory
2198   // File name will be specified by : sector, side, runid and eventnumber
2199   
2200   if(GetEventProcessed())
2201   {
2202     AliInfo("Write histos to file");
2203     Char_t name[256];
2204     sprintf(name,"SIDE_%i_SECTOR_%02i_RUN_%05i_EventID_%06i.root",(GetLastSector()/18),(GetLastSector()%18),fRunId,fEventNumber);
2205     TFile* f = new TFile(name,"recreate");
2206     for(Int_t i =0; i<fHistList->GetEntries(); i++)
2207     {
2208       if(((TH1*)fHistList->At(i))!=0)
2209       {
2210         ((TH1*)fHistList->At(i))->Write();
2211       }
2212     }
2213     f->ls();
2214     f->Close();
2215   }
2216   else
2217   {
2218     AliError("No Event Processed : Chose Format , File and push 'Next Event' ");
2219   }
2220 }
2221
2222
2223 //__________________________________________________________________
2224 TH1* AliTPCMonitor::GetHisto(char* histname) 
2225 {
2226   
2227   // Returns histogram specified by histname
2228   // check available names for histos in CreateHistos()
2229   
2230   TH1* hist = 0;
2231   if((TH1*)fHistList->FindObject(histname))
2232   {
2233     hist = (TH1*)fHistList->FindObject(histname);
2234   }
2235   else
2236   {
2237     cout << " AliTPCMonitor::GetHisto :: Can not find histo with name " << histname << endl;
2238   }
2239   return hist ;
2240 }