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