]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSOnlineMonitor.cxx
Raw Sdigits energy converted to GeV.
[u/mrichter/AliRoot.git] / PHOS / AliPHOSOnlineMonitor.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 //_________________________________________________________________________
19 // Class intended to perform online monitoring of PHOS beamtests
20 // Being constructed, produces menu with list of available histograms to fill
21 // Once histograms are selected, button "Go" should be pressed to start scan of data.
22 // Prepared histograms will be periodically updated during scan of the data.
23 // Note:
24 // 1. To plot most of histograms, a "Connection Table", relating ADC signal index and AbsId 
25 //    of PHOS crystal, should be created beforehand. To do this, call macro
26 //       $ALICE_ROOT/PHOS/macros/BeamTest/MakeConTableDB.C 
27 //    with apropriate number of raws and columns of prototype.
28 // 2. To perform reconstruction and e.g. invariant mass analysis, a "Calibration Database"
29 //    should be created beforehand. To do this, call macro 
30 //       $ALICE_ROOT/PHOS/macros/BeamTest/MakeConTableDB.C 
31 //    to read calibration parameters from file or use AliPHOSCalibrator to calculate 
32 //    pedestal and gains.
33 // 3. Once histograms are filled with "Go" method, they can be written to file
34 //    with WriteHisto("Filename.root") method. 
35 //
36 //*-- Author : D.Peressounko (RRC KI) after A.V. Kuryakin, (Sarov) 
37 //////////////////////////////////////////////////////////////////////////////
38
39 // --- ROOT system ---
40 #include "TROOT.h"
41 #include "TClonesArray.h"
42 #include "TCanvas.h"
43 #include "TH1.h"
44 #include "TH2.h"
45 #include "TGroupButton.h"
46 #include "TFile.h"
47 #include "TSystem.h"
48
49
50 // --- Standard library ---
51 #include "TBenchmark.h"
52 #include "Riostream.h"
53
54 // --- AliRoot header files ---
55 #include "AliPHOSOnlineMonitor.h"
56 #include "AliPHOSConTableDB.h"
57 #include "AliPHOSGeometry.h"
58 #include "AliRawReaderDateV3.h"
59 #include "AliRawEventHeaderBase.h"
60 #include "AliPHOSRawStream2004.h"
61 #include "AliPHOSDigit.h"
62 #include "AliPHOSGetterLight.h"  
63 #include "AliPHOSClusterizerv1.h"  
64 #include "AliPHOSTrackSegmentMakerv1.h"  
65 #include "AliPHOSPIDv1.h"  
66 #include "AliPHOSCalibrManager.h" 
67 #include "AliPHOSCalibrationDB.h"
68
69 ClassImp(AliPHOSOnlineMonitor)
70
71
72 //____________________________________________________________________________ 
73 AliPHOSOnlineMonitor::AliPHOSOnlineMonitor(): 
74   TDialogCanvas("PHOS","PHOS",150,300),
75   fScanPed(kFALSE),
76   fScanSig(kFALSE),
77   fReconstruct(kFALSE),
78   fNevents(0),
79   fNUpdate(1000),
80   fCanvasList(new TList),
81   fHistosList(new TList),
82   fInputFile(),
83   fGeom(0),
84   fcdb(0)
85 {
86   MakeButtons() ;
87   Modified(kTRUE);
88   Update();
89   SetEditable(kFALSE);
90   
91   //add this TFitPanel to the list of cleanups such that in case
92   //the referenced object is deleted, its pointer be reset
93   gROOT->GetListOfCleanups()->Add(this);
94   
95   fRefObject = this;
96   fRefPad    = (TPad*)gROOT->GetSelectedPad();
97
98   fGeom = AliPHOSGeometry::GetInstance("IHEP","") ;
99 }
100 //____________________________________________________________________________ 
101 AliPHOSOnlineMonitor::AliPHOSOnlineMonitor(const char * inputfile): 
102   TDialogCanvas("PHOS","PHOS",150,300),
103   fScanPed(kFALSE),
104   fScanSig(kFALSE),
105   fReconstruct(kFALSE),
106   fNevents(0),
107   fNUpdate(1000),
108   fCanvasList(new TList),
109   fHistosList(new TList),
110   fInputFile(inputfile),
111   fGeom(0),
112   fcdb(0)
113 {
114   MakeButtons() ;
115   Modified(kTRUE);
116   Update();
117   SetEditable(kFALSE);
118   
119   //add this TFitPanel to the list of cleanups such that in case
120   //the referenced object is deleted, its pointer be reset
121   gROOT->GetListOfCleanups()->Add(this);
122   
123   fRefObject = this;
124   fRefPad    = (TPad*)gROOT->GetSelectedPad();
125
126   fGeom = AliPHOSGeometry::GetInstance("IHEP","") ;
127 }
128
129 //____________________________________________________________________________ 
130 AliPHOSOnlineMonitor::AliPHOSOnlineMonitor(const AliPHOSOnlineMonitor & /*rhs*/):
131   TDialogCanvas(),
132   fScanPed(kFALSE),
133   fScanSig(kFALSE),
134   fReconstruct(kFALSE),
135   fNevents(0),
136   fNUpdate(0),
137   fCanvasList(0),
138   fHistosList(0),
139   fInputFile(),
140   fGeom(0),
141   fcdb(0)
142 {
143   Fatal("AliPHOSOnlineMonitor", "not implemented");
144 }
145
146 //____________________________________________________________________________ 
147 AliPHOSOnlineMonitor & AliPHOSOnlineMonitor::operator = (const AliPHOSOnlineMonitor &)
148 {
149   Fatal("operator = ", "not implemented");
150   return *this;
151 }
152
153 //____________________________________________________________________________ 
154 AliPHOSOnlineMonitor::~AliPHOSOnlineMonitor()
155 {
156   //Obvious, but unevoidable comment for destructor: cleans up everething.
157   TIter nextCanvas(fCanvasList);
158   TCanvas * c ;
159   while((c=(TCanvas*)nextCanvas()))
160     delete c ;
161   delete fCanvasList ;
162   
163   TIter nextHisto(fHistosList);
164   TH1D * h ;
165   while((h=(TH1D*)nextHisto()))
166     delete h ;
167   delete fHistosList ; 
168
169   if(fcdb)
170     delete fcdb ;
171
172 }
173 //____________________________________________________________________________ 
174 void AliPHOSOnlineMonitor::MakeButtons(void){
175   //Make buttons on graphical menu
176   Int_t nButtons = 16;
177   TGroupButton * b ;
178   Float_t xmin = 0.0;
179   Float_t ymin = 0.01;
180   Float_t xmax = 0.99;
181   Float_t ymax = 0.99;
182   Float_t dy = (ymax-ymin)/nButtons ;
183
184   Float_t y2=ymax ;
185   Float_t y1=y2-dy ;
186   b = new TGroupButton("APPLY","Triggers","",xmin,y1,xmax,y2);
187   b->Draw();
188   y2=y1 ;
189   y1=y1-dy ;
190   b = new TGroupButton("APPLY","Pedestals","",xmin,y1,xmax,y2);
191   b->Draw();
192   y2=y1 ;
193   y1=y1-dy ;
194   b = new TGroupButton("APPLY","Spectrum all","",xmin,y1,xmax,y2);
195   b->Draw();
196   y2=y1 ;
197   y1=y1-dy ;
198   b = new TGroupButton("APPLY","Spectrum g","",xmin,y1,xmax,y2);
199   b->Draw();
200   y2=y1 ;
201   y1=y1-dy ;
202   b = new TGroupButton("APPLY","Inv Mass","",xmin,y1,xmax,y2);
203   b->Draw();
204   for(Int_t i=1; i<=5; i++){
205     y2=y1 ;
206     y1=y1-dy ;
207     char name[10] ;
208     sprintf(name,"Edep(ADC) %d",i) ;
209     b = new TGroupButton("APPLY",name,"",xmin,y1,xmax,y2);
210     b->Draw();
211   }
212   for(Int_t i=1; i<=5; i++){
213     y2=y1 ;
214     y1=y1-dy ;
215     char name[10] ;
216     sprintf(name,"Edep(Cal) %d",i) ;
217     b = new TGroupButton("APPLY",name,"",xmin,y1,xmax,y2);
218     b->Draw();
219   }
220   y2=y1 ;
221   y1=y1-dy ;
222   b = new TGroupButton("APPLY","Go","",xmin,y1,xmax,y2);
223   b->SetTextColor(2);
224   b->Draw();
225 }
226 //____________________________________________________________________________ 
227 void AliPHOSOnlineMonitor::Apply(const char *action){
228   //Function to handle button actions
229
230   TDialogCanvas::Apply() ;
231
232   TObject *obj;
233   TGroupButton *button;
234   TIter next(fPrimitives);
235   
236   if (!strcmp(action,"Triggers")) {
237     DrawTriggers() ;
238   }
239
240   if (!strcmp(action,"Pedestals")) {
241     DrawPedestals() ;
242   }
243   if (!strcmp(action,"Spectrum all")) {
244     DrawSpectrum("all") ;
245   }
246   if (!strcmp(action,"Spectrum g")) {
247     DrawSpectrum("gamma") ;
248   }
249   if (!strcmp(action,"Inv Mass")) {
250     DrawMinv() ;
251   }
252   if(strstr(action,"Edep")){
253     Int_t n ;
254     char tmp[10] ;
255     sscanf(action,"%s %d",tmp,&n) ;
256     char opt[5]="" ;
257     if(strstr(action,"Cal"))
258       sprintf(opt,"Cal") ;
259     DrawEdep(n,opt) ;
260   }
261   if (!strcmp(action,"Go")) {
262     Go() ;
263   }
264
265   //Mark button as pressed
266   if(strcmp(action,"Go")){ //Do not mark "Go" button
267     while ((obj = next())) {
268       if (obj->InheritsFrom(TGroupButton::Class())) {
269         button = (TGroupButton*)obj;
270         if(!strcmp(button->GetTitle(),action)){
271           if (button->GetBorderMode() > 0){
272             button->SetBorderMode(-1) ;
273             button->Modified(kTRUE);
274           }
275         }
276       }
277     }
278   }
279 }
280 //____________________________________________________________________________ 
281 void  AliPHOSOnlineMonitor::SetInputFile(const char * filename){
282   //close previously opened
283   
284   fInputFile = filename ;
285 }
286 //____________________________________________________________________________ 
287 void  AliPHOSOnlineMonitor::DrawPedestals(){
288   //Prepare canvas and histograms for drawing pedestals
289
290   TIter nextCanvas(fCanvasList);
291   TCanvas * c ;
292   Bool_t exists = kFALSE ;
293   while((c=(TCanvas*)nextCanvas())){
294     if(!strcmp(c->GetName(),"Pedestals")){
295       exists = kTRUE ;
296       break;
297     }
298   }
299   if(!exists){
300     c = new TDialogCanvas("Pedestals","Pedestals",300,200) ;
301     fCanvasList->AddLast(c) ;
302   }
303   
304   TIter nextHisto(fHistosList);
305   TH1D * h ;
306   exists = kFALSE ;
307   while((h=(TH1D*)nextHisto())){
308     if(!strcmp(h->GetName(),"hPedestals")){
309       exists = kTRUE ;
310       break;
311     }
312   }
313   if(!exists){
314     h = new TH1D("hPedestals","Pedestals per event",fGeom->GetNModules()*fGeom->GetNCristalsInModule(),0.,
315                  1.*fGeom->GetNModules()*fGeom->GetNCristalsInModule()) ;
316     fHistosList->AddLast(h) ;
317   }
318   
319   c->cd() ;
320   h->SetStats(0) ;
321   h->Draw() ;
322   fScanPed = kTRUE ; //We will scan pedestals
323
324 }
325 //____________________________________________________________________________ 
326 void  AliPHOSOnlineMonitor::DrawTriggers(){
327   //Prepare canvas and histogram for drawing triggers
328
329   TIter nextCanvas(fCanvasList);
330   TCanvas * c ;
331   Bool_t exists = kFALSE ;
332   while((c=(TCanvas*)nextCanvas())){
333     if(!strcmp(c->GetName(),"Triggers")){
334       exists = kTRUE ;
335       break;
336     }
337   }
338   if(!exists){
339     c = new TDialogCanvas("Triggers","Triggers",200,200) ;
340     fCanvasList->AddLast(c) ;
341   }
342   
343   TIter nextHisto(fHistosList);
344   TH1D * h ;
345   exists = kFALSE ;
346   while((h=(TH1D*)nextHisto())){
347     if(!strcmp(h->GetName(),"hTriggers")){
348       exists = kTRUE ;
349       break;
350     }
351   }
352   if(!exists){
353     h = new TH1D("hTriggers","Triggers",2,0.,2.) ;
354     fHistosList->AddLast(h) ;
355   }
356   //Make Labels
357   h->SetBit(TH1::kCanRebin);  
358   h->Fill("LED",0.0000001) ; 
359   h->Fill("PUL",0.0000001) ; 
360   h->Fill("PED",0.0000001) ; 
361   h->Fill("NEL",0.0000001) ; 
362   h->Fill("WEL",0.0000001) ; 
363   h->Fill("SOB",0.0000001) ; 
364   h->Fill("EOB",0.0000001) ; 
365   h->Fill("wrong",0.0000001) ;
366   h->LabelsOption("h");
367   h->LabelsDeflate();
368   h->SetStats(0) ;
369   c->cd() ;
370   h->Draw() ;
371 }
372 //____________________________________________________________________________ 
373 void  AliPHOSOnlineMonitor::DrawSpectrum(const char * opt){
374   //Prepare canvas and histograms for drawing spectra of all reconstructed particles or photons
375
376   TString name("Spectrum") ;
377   name+=opt ;
378   
379   TIter nextCanvas(fCanvasList);
380   TCanvas * c ;
381   Bool_t exists = kFALSE ;
382   while((c=(TCanvas*)nextCanvas())){
383     if(!strcmp(c->GetName(),name.Data())){
384       exists = kTRUE ;
385       break;
386     }
387   }
388   if(!exists){
389     c = new TDialogCanvas(name,name,250,300) ;
390     fCanvasList->AddLast(c) ;
391   }
392   
393   TIter nextHisto(fHistosList);
394   TH1D * h ;
395   exists = kFALSE ;
396   name.Prepend("h") ;
397   while((h=(TH1D*)nextHisto())){
398     if(!strcmp(h->GetName(),name.Data())){
399       exists = kTRUE ;
400       break;
401     }
402   }
403   if(!exists){
404     h = new TH1D(name,name,100,0.,100.) ;
405     fHistosList->AddLast(h) ;
406   }
407   
408   h->SetStats(0) ;
409   c->cd() ;
410   h->Draw() ;
411   fReconstruct = kTRUE ;
412   fScanSig = kTRUE ; //We will scan pedestals
413
414 }
415 //____________________________________________________________________________ 
416 void  AliPHOSOnlineMonitor::DrawMinv(){
417   TIter nextCanvas(fCanvasList);
418   TCanvas * c ;
419   Bool_t exists = kFALSE ;
420   while((c=(TCanvas*)nextCanvas())){
421     if(!strcmp(c->GetName(),"InvMass")){
422       exists = kTRUE ;
423       break;
424     }
425   }
426   if(!exists){
427     c = new TDialogCanvas("InvMass","Invariant mass",300,200) ;
428     fCanvasList->AddLast(c) ;
429   }
430   
431   TIter nextHisto(fHistosList);
432   TH1D * h ;
433   exists = kFALSE ;
434   while((h=(TH1D*)nextHisto())){
435     if(!strcmp(h->GetName(),"hInvMass")){
436       exists = kTRUE ;
437       break;
438     }
439   }
440   if(!exists){
441     h = new TH1D("hInvMass","hInvMass",1000,0.,1.0) ;
442     fHistosList->AddLast(h) ;
443   }
444   
445   c->cd() ;
446   h->Draw() ;
447   h->SetStats(0) ;
448   fReconstruct = kTRUE ;
449   fScanSig = kTRUE ; //We will scan pedestals
450 }
451 //____________________________________________________________________________ 
452 void  AliPHOSOnlineMonitor::DrawEdep(Int_t mod,const char * opt){
453   char name[15] ;
454   sprintf(name,"Edep%s %d",opt,mod) ;
455
456   TIter nextCanvas(fCanvasList);
457   TCanvas * c ;
458   Bool_t exists = kFALSE ;
459   while((c=(TCanvas*)nextCanvas())){
460     if(!strcmp(c->GetName(),name)){
461       exists = kTRUE ;
462       break;
463     }
464   }
465   if(!exists){
466     c = new TDialogCanvas(name,name,300,200) ;
467     fCanvasList->AddLast(c) ;
468   }
469   
470   TIter nextHisto(fHistosList);
471   TH2D * h ;
472   exists = kFALSE ;
473   sprintf(name,"hEdep%s%d",opt,mod) ;
474   while((h=(TH2D*)nextHisto())){
475     if(!strcmp(h->GetName(),name)){
476       exists = kTRUE ;
477       break;
478     }
479   }
480   if(!exists){
481     h = new TH2D(name,name,fGeom->GetNPhi(),0.,1.*fGeom->GetNPhi(),fGeom->GetNZ(),0.,1.*fGeom->GetNZ()) ;
482     fHistosList->AddLast(h) ;
483   }
484   
485   c->cd() ;
486   h->Draw("col") ;
487   h->SetStats(0) ;
488   fScanSig = kTRUE ; //We will scan signal events
489   if(strstr(opt,"Cal"))
490     fReconstruct = kTRUE ;
491  
492 }
493 //____________________________________________________________________________ 
494 void  AliPHOSOnlineMonitor::ScanPedestals(TClonesArray * digits){
495   //This method is called for events with PED trigger
496   //We fill bins with ADC values
497
498   TH1D * h = (TH1D*)gROOT->FindObjectAny("hPedestals");
499   if(!h){
500     Error("ScanPedestals","Can not fild histogram hPedestals") ;
501     return ;
502   }
503   for(Int_t i=0; i<digits->GetEntriesFast(); i++){
504     AliPHOSDigit * dig = static_cast<AliPHOSDigit*>(digits->At(i)) ;
505     h->AddBinContent(dig->GetId(),dig->GetAmp()) ;
506   }
507 }
508 //____________________________________________________________________________ 
509 void  AliPHOSOnlineMonitor::ScanEdep(TClonesArray * digits){
510   //Fill 2D distribution of ADC values in NEL and WEL events
511   AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
512   AliPHOSCalibrationDB *cdb = 0 ;
513   if(gime)
514     cdb = gime->CalibrationDB() ;
515   Int_t mod = 0 ;
516   char name[15] ;
517   TH2D * h = 0 ;
518   TH2D * hCal = 0 ;
519   for(Int_t i=0; i<digits->GetEntriesFast(); i++){
520     AliPHOSDigit * dig = static_cast<AliPHOSDigit*>(digits->At(i)) ;
521     Int_t relId[4] ;
522     fGeom->AbsToRelNumbering(dig->GetId(),relId) ;
523     if(mod != relId[0]){ //new module, look for histograms
524       mod = relId[0] ;
525       sprintf(name,"hEdep%d",mod) ;
526       h = (TH2D*)gROOT->FindObjectAny(name);
527       sprintf(name,"hEdepCal%d",mod) ;
528       hCal = (TH2D*)gROOT->FindObjectAny(name);
529     }
530     if(h)
531       h->Fill(relId[2]-0.1,relId[3]-0.1,1.*dig->GetAmp()) ;
532     if(hCal)
533       hCal->Fill(relId[2]-0.1,relId[3]-0.1,cdb->Calibrate(dig->GetAmp(),dig->GetId())) ;
534   }
535 }
536 //____________________________________________________________________________ 
537 void  AliPHOSOnlineMonitor::ScanRecon(TClonesArray * recParticles){
538   if(!recParticles || recParticles->GetEntries()==0) return ;
539
540   TH1D* hSpectr   = (TH1D*)gROOT->FindObjectAny("hSpectrumall");
541   TH1D* hSpectrGam= (TH1D*)gROOT->FindObjectAny("hSpectrumgamma");
542   TH1D* hInvMass  = (TH1D*)gROOT->FindObjectAny("hInvMass");
543   for(Int_t i=0; i<recParticles->GetEntriesFast() ; i++){
544     AliPHOSRecParticle * p = (AliPHOSRecParticle *)recParticles->At(i) ;
545     if(hSpectr)hSpectr->Fill(p->Energy()) ;
546     if(hSpectrGam && p->IsPhoton())hSpectrGam->Fill(p->Energy()) ;
547     if(hInvMass){
548       for(Int_t j=i+1; j<recParticles->GetEntriesFast() ; j++){
549         AliPHOSRecParticle * p2 = (AliPHOSRecParticle *)recParticles->At(j) ;
550         Double_t e = p->Energy() + p2->Energy() ;
551         Double_t x = p->Px() + p2->Px() ;
552         Double_t y = p->Py() + p2->Py() ;
553         Double_t z = p->Pz() + p2->Pz() ;
554         Double_t m = e*e-x*x-y*y-z*z ;
555         hInvMass->Fill(m>0?TMath::Sqrt(m): 0. ) ;
556       }
557     }
558   }
559 }
560 //____________________________________________________________________________ 
561 void  AliPHOSOnlineMonitor::ScanTrigger(Int_t trig){
562   //Fills trigger distribution
563
564   TH1D * h = (TH1D*)gROOT->FindObjectAny("hTriggers");
565   if(!h) return ;
566   switch(trig){
567   case AliPHOSRawStream2004::kLED : h->Fill("LED",1.) ; break ;
568   case AliPHOSRawStream2004::kPUL : h->Fill("PUL",1.) ; break ;
569   case AliPHOSRawStream2004::kPED : h->Fill("PED",1.) ; break ;
570   case AliPHOSRawStream2004::kNEL : h->Fill("NEL",1.) ; break ;
571   case AliPHOSRawStream2004::kWEL : h->Fill("WEL",1.) ; break ;
572   case AliPHOSRawStream2004::kSOB : h->Fill("SOB",1.) ; break ;
573   case AliPHOSRawStream2004::kEOB : h->Fill("EOB",1.) ; break ;
574   default : h->Fill("wrong",1.) ;
575   }
576 }
577 //____________________________________________________________________________ 
578 void  AliPHOSOnlineMonitor::SetConTableDB(const char * filename){
579   //Read ConnectionTableDB from file
580   TFile * file = new TFile(filename) ;
581   AliPHOSConTableDB * tmp = (AliPHOSConTableDB*)file->Get("AliPHOSConTableDB") ;
582   fcdb = new AliPHOSConTableDB(*tmp) ;
583   file->Close() ; 
584 }
585 //____________________________________________________________________________ 
586 void  AliPHOSOnlineMonitor::Go(){
587   //Perform scan of curent event
588   gBenchmark->Start("PHOSOnlineMon"); 
589
590   //First test if we need "Connection table" then open it
591   if(!fcdb){
592     SetConTableDB() ;
593     if(fcdb){
594       Info("Go","Read Connection table from file \"ConTableDB.root\"") ;
595     }else{
596       Error("Go","Please, set connection table with SetConTableDB() method") ;
597       return ;
598     }
599   }
600
601   AliPHOSGetterLight * gime = AliPHOSGetterLight::Instance("PHOS","On Flight") ;
602
603   //Configure CalibrManager to read data from file
604   //Create calibration database and read it
605   AliPHOSCalibrationDB * calibDB = 0 ;
606   if(fScanSig || fReconstruct){ //We will ned calibration parameters
607     AliPHOSCalibrManager::GetInstance("CalibrDB.root","root") ;
608     //If we configured manager to read from ASCII file, 
609     //give him connection table. OK, it will not harm in any case.
610     AliPHOSCalibrManager::GetInstance()->SetConTable(fcdb) ;
611     
612     calibDB = new AliPHOSCalibrationDB("OnLine") ;
613     calibDB->GetParameters() ; //Read parameters using Manager
614     gime->SetCalibrationDB(calibDB) ;
615   }
616   
617   //Now open data file
618   AliRawReaderDateV3 *rawReader = new AliRawReaderDateV3(fInputFile) ; 
619   rawReader->RequireHeader(kFALSE);
620   AliPHOSRawStream2004     *rawStream = new AliPHOSRawStream2004(rawReader) ;
621   rawStream->SetConTableDB(fcdb) ;
622   
623   TClonesArray * digits = gime->Digits() ;
624   TClonesArray * recParticles = gime->RecParticles() ;
625   AliPHOSClusterizerv1* clu = 0 ;
626   AliPHOSTrackSegmentMakerv1 * tsm = 0 ;
627   AliPHOSPIDv1 * pid = 0 ;
628   if(fReconstruct){ //We will need calibation parameters
629     clu = new AliPHOSClusterizerv1("PHOS","On Flight") ;
630     clu->SetWriting(0) ; //Do not write to file
631     clu->SetEmcMinE(0.05) ;  //Minimal energy of the digit
632     clu->SetEmcLocalMaxCut(0.05) ; //Height of local maximum over environment
633     clu->SetEmcClusteringThreshold(0.2) ; //Minimal energy to start cluster
634 //    clu->SetUnfolding(kFALSE) ; //Do not unfold
635     tsm = new AliPHOSTrackSegmentMakerv1("PHOS","On Flight") ;
636     tsm->SetWriting(0) ; //Do not write to file
637     pid = new AliPHOSPIDv1("PHOS","On Flight") ;
638     pid->SetWriting(0) ; //Do not write to file    
639   }
640   
641   fNevents=0 ;
642   //Scan all event in file
643   printf("      ") ;
644   while(rawReader->NextEvent()){
645     //Is it PHYSICAL event
646     if(rawReader->GetType() == AliRawEventHeaderBase::kPhysicsEvent){
647       fNevents++ ;
648       if(fNevents%100 ==0){
649         printf("\b\b\b\b\b\b%6d",fNevents) ;
650       }
651       if(rawStream->ReadDigits(digits)){
652         
653         //Test trigger
654         //Pedestal Event
655         ScanTrigger(rawStream->GetTrigger()) ;
656         if(rawStream->IsPEDevent() && fScanPed){
657           ScanPedestals(digits) ;
658         }
659         if((rawStream->IsNELevent() || rawStream->IsWELevent()) && fScanSig){
660           ScanEdep(digits) ;
661           if(fReconstruct){
662             gime->Clusterizer()->Exec("") ;
663             gime->TrackSegmentMaker()->Exec("") ;
664             gime->PID()->Exec("") ;
665             ScanRecon(recParticles) ;
666           }
667         }
668       }
669       
670       if(fNevents%fNUpdate == 0 ){ //upate all histograms       
671         TIter nextCanvas(fCanvasList);
672         TCanvas * c ;
673         while((c=(TCanvas*)nextCanvas())){
674           c->Modified() ;
675           c->Update() ;
676         }
677       }
678       gSystem->ProcessEvents(); 
679     }
680     //    if(fNevents>=200)break ;
681   }
682   printf("\n") ;
683   gBenchmark->Stop("PHOSOnlineMon");
684   Float_t time = gBenchmark->GetCpuTime("PHOSOnlineMon") ;
685   printf("took %f seconds for scanning, i.e. %f seconds per event %d  \n",
686          time,time/fNevents,fNevents) ; 
687   
688   //Update canvas with histograms at the end
689   TIter nextCanvas(fCanvasList);
690   TCanvas * c ;
691   while((c=(TCanvas*)nextCanvas())){
692     c->Modified(kTRUE) ;
693   }
694   
695   if(clu)delete clu ;
696   if(tsm)delete tsm ;
697   if(pid)delete pid ;
698   printf("delete 1 \n") ;
699   if(calibDB) delete calibDB ;
700   delete rawStream ;
701   delete rawReader ;
702 }
703 //____________________________________________________________________________ 
704 void  AliPHOSOnlineMonitor::Clean(){
705   //Cleans content of all histograms
706
707   TIter nextHisto(fHistosList);
708   TH1D * h ;
709   while((h=(TH1D*)nextHisto())){
710     h->Reset("ISE") ;
711   }
712   TIter nextCanvas(fCanvasList);
713   TCanvas * c ;
714   while((c=(TCanvas*)nextCanvas())){
715     c->Modified() ;
716   }
717 }
718 //____________________________________________________________________________ 
719 void  AliPHOSOnlineMonitor::Reset(){
720   //delets all canvas and histograms,
721   //marks buttons as unpressed
722
723   TIter nextHisto(fHistosList);
724   TH1D * h ;
725   while((h=(TH1D*)nextHisto())){
726     fHistosList->Remove(h) ;
727     delete h ;
728   }
729   TIter nextCanvas(fCanvasList);
730   TCanvas * c ;
731   while((c=(TCanvas*)nextCanvas())){
732     fCanvasList->Remove(c) ;
733     delete c ;
734   }
735   TObject *obj;
736   TGroupButton *button;
737   TIter next(fPrimitives);
738   
739   //Mark buttons as anpressed
740   while ((obj = next())) {
741     if (obj->InheritsFrom(TGroupButton::Class())) {
742       button = (TGroupButton*)obj;
743       if (button->GetBorderMode() < 0){
744         button->SetBorderMode(1) ;
745         button->Modified(kTRUE);
746       }
747     }
748   }
749   
750 }
751 //____________________________________________________________________________ 
752 void  AliPHOSOnlineMonitor::WriteHistograms(const char * filename){
753   //Write filled histograms to file
754   TFile * file = new TFile(filename,"Update") ;
755   file->cd() ;
756   TIter nextHisto(fHistosList);
757   TH1 * h ;
758   while((h=(TH1*)nextHisto())){
759     h->Write(0,TObject::kOverwrite) ;
760   }
761   file->Close() ; 
762 }