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