]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MONITOR/AliMonitorProcess.cxx
monitoring of HLT hough transformation added
[u/mrichter/AliRoot.git] / MONITOR / AliMonitorProcess.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 ///////////////////////////////////////////////////////////////////////////////
19 //                                                                           //
20 //  This is the class for perfoming the monitoring process.                  //
21 //  It checks if a raw data file exists, loops over the events in the raw    //
22 //  data file, reconstructs TPC and ITS clusters and tracks, fills the       //
23 //  monitor histograms and sends the updated histograms to the clients.      //
24 //  Then the raw data file is deleted and it waits for a new file.           //
25 //                                                                           //
26 ///////////////////////////////////////////////////////////////////////////////
27
28
29 #include "AliMonitorProcess.h"
30 #include "AliMonitorTPC.h"
31 #include "AliMonitorITS.h"
32 #include "AliMonitorV0s.h"
33 #include "AliMonitorHLT.h"
34 #include "AliMonitorHLTHough.h"
35 #include "AliRawReaderRoot.h"
36 #include "AliLoader.h"
37 #include "AliRun.h"
38 #include "AliTPC.h"
39 #include "AliTPCclustererMI.h"
40 #include "AliTPCtrackerMI.h"
41 #include "AliITS.h"
42 #include "AliITSclustererV2.h"
43 #include "AliITStrackerV2.h"
44 #include "AliITSLoader.h"
45 #include "AliV0vertexer.h"
46 #include <TGrid.h>
47 #include <TSocket.h>
48 #include <TServerSocket.h>
49 #include <TMessage.h>
50 #include <TGridResult.h>
51 #include <TROOT.h>
52 #ifdef ALI_HLT
53 #include <AliLevel3.h>
54 #include <AliL3Transform.h>
55 #include <AliL3StandardIncludes.h>
56 #include <AliL3HoughMaxFinder.h>
57 #include <AliL3Hough.h>
58 #endif
59
60 ClassImp(AliMonitorProcess) 
61
62
63 const Int_t AliMonitorProcess::fgkPort = 9327;
64
65
66 //_____________________________________________________________________________
67 AliMonitorProcess::AliMonitorProcess(const char* alienDir,
68                                      const char* fileNameGalice)
69 {
70 // initialize the monitoring process and the monitor histograms
71
72   fGrid = TGrid::Connect("alien", gSystem->Getenv("USER"));
73   if (!fGrid || fGrid->IsZombie() || !fGrid->IsConnected()) {
74     delete fGrid;
75     Fatal("AliMonitorProcess", "could not connect to alien");
76   }
77   fGrid->cd(alienDir);
78   fLogicalFileName = "";
79   fFileName = "";
80
81   fRunLoader = AliRunLoader::Open(fileNameGalice);
82   if (!fRunLoader) Fatal("AliMonitorProcess", 
83                          "could not get run loader from file %s", 
84                          fileNameGalice);
85
86   fRunLoader->CdGAFile();
87   fTPCParam = AliTPC::LoadTPCParam(gFile);
88   if (!fTPCParam) Fatal("AliMonitorProcess", "could not load TPC parameters");
89
90   fRunLoader->LoadgAlice();
91   gAlice = fRunLoader->GetAliRun();
92   if (!gAlice) Fatal("AliMonitorProcess", "no gAlice object found");
93   AliITS* its = (AliITS*) gAlice->GetModule("ITS");
94   if (!its) Fatal("AliMonitorProcess", "no ITS detector found");
95   fITSgeom = its->GetITSgeom();
96   if (!fITSgeom) Fatal("AliMonitorProcess", "could not load ITS geometry");
97
98 #ifdef ALI_HLT
99 // Init TPC parameters for HLT
100   Bool_t isinit=AliL3Transform::Init(const_cast<char*>(fileNameGalice),kTRUE);
101   if(!isinit){
102     cerr << "Could not create transform settings, please check log for error messages!" << endl;
103     return;
104   }
105 #endif
106
107   fRunNumber = 0;
108   fSubRunNumber = 0;
109   fNEvents = 0;
110   fNEventsMin = 2;
111   fWriteHistoList = kFALSE;
112
113   fTopFolder = new TFolder("Monitor", "monitor histograms");
114   fTopFolder->SetOwner(kTRUE);
115
116   fMonitors.Add(new AliMonitorTPC(fTPCParam));
117   fMonitors.Add(new AliMonitorITS(fITSgeom));
118   fMonitors.Add(new AliMonitorV0s);
119 #ifdef ALI_HLT
120   fMonitors.Add(new AliMonitorHLT(fTPCParam));
121   fMonitors.Add(new AliMonitorHLTHough(fTPCParam));
122 #endif
123
124   for (Int_t iMonitor = 0; iMonitor < fMonitors.GetEntriesFast(); iMonitor++) {
125     ((AliMonitor*) fMonitors[iMonitor])->CreateHistos(fTopFolder);
126   }
127
128   TIterator* iFolder = fTopFolder->GetListOfFolders()->MakeIterator();
129   while (TFolder* folder = (TFolder*) iFolder->Next()) folder->SetOwner(kTRUE);
130   delete iFolder;
131
132   fFile = TFile::Open("monitor_tree.root", "RECREATE");
133   if (!fFile || !fFile->IsOpen()) {
134     Fatal("AliMonitorProcess", "could not open file for tree");
135   }
136   fTree = new TTree("MonitorTree", "tree for monitoring");
137   for (Int_t iMonitor = 0; iMonitor < fMonitors.GetEntriesFast(); iMonitor++) {
138     ((AliMonitor*) fMonitors[iMonitor])->CreateBranches(fTree);
139   }
140   gROOT->cd();
141
142   fServerSocket = new TServerSocket(fgkPort, kTRUE);
143   fServerSocket->SetOption(kNoBlock, 1);
144   fDisplaySocket = NULL;
145   CheckForConnections();
146 #ifdef ALI_HLT
147   fHLT = NULL;
148 #endif
149
150   SetStatus(kStopped);
151   fStopping = kFALSE;
152 }
153
154 //_____________________________________________________________________________
155 AliMonitorProcess::AliMonitorProcess(const AliMonitorProcess& process) :
156   TObject(process)
157 {
158   Fatal("AliMonitorProcess", "copy constructor not implemented");
159 }
160
161 //_____________________________________________________________________________
162 AliMonitorProcess& AliMonitorProcess::operator = (const AliMonitorProcess& 
163                                                   /*process*/)
164 {
165   Fatal("operator =", "assignment operator not implemented");
166   return *this;
167 }
168
169 //_____________________________________________________________________________
170 AliMonitorProcess::~AliMonitorProcess()
171 {
172 // clean up
173
174   fMonitors.Delete();
175   delete fTopFolder;
176   delete fRunLoader;
177
178   delete fServerSocket;
179   fSockets.Delete();
180   delete fDisplaySocket;
181
182   fGrid->Close();
183   delete fGrid;
184
185   fFile->Close();
186   delete fFile;
187   gSystem->Unlink("monitor_tree.root");
188
189 #ifdef ALI_HLT
190   delete fHLT;
191 #endif
192 }
193
194
195 //_____________________________________________________________________________
196 const char* AliMonitorProcess::GetRevision()
197 {
198   return "$Revision$";
199 }
200
201
202 //_____________________________________________________________________________
203 void AliMonitorProcess::SetStatus(EStatus status)
204 {
205 // set the current status and process system events
206
207   fStatus = status;
208   gSystem->ProcessEvents();
209 }
210
211
212 //_____________________________________________________________________________
213 void AliMonitorProcess::Run()
214 {
215 // run the monitor process: 
216 //  check for a raw data file, process the raw data file and delete it
217
218   fStopping = kFALSE;
219
220   while (!fStopping) {
221     SetStatus(kWaiting);
222     while (!CheckForNewFile()) {
223       CheckForConnections();
224       SetStatus(kWaiting);
225       if (fStopping) break;
226       gSystem->Sleep(10);
227     }
228     if (fStopping) break;
229
230     ProcessFile();
231   }
232
233   WriteHistos();
234
235   fStopping = kFALSE;
236   SetStatus(kStopped);
237 }
238
239
240 //_____________________________________________________________________________
241 void AliMonitorProcess::Stop()
242 {
243 // set the fStopping flag to terminate the monitor process after the current
244 // event was processed
245
246   if (GetStatus() != kStopped) fStopping = kTRUE;
247 }
248
249
250 //_____________________________________________________________________________
251 void AliMonitorProcess::ProcessFile(const char* fileName)
252 {
253 // create a file with monitor histograms for a single file
254
255   if (GetStatus() != kStopped) {
256     Error("ProcessFile", "ProcessFile can not be called"
257           " while the monitor process is running");
258     return;
259   }
260
261   fFileName = fileName;
262   Int_t nEventMin = fNEventsMin;
263   fNEventsMin = 1;
264   ProcessFile();
265   WriteHistos();
266   fNEventsMin = nEventMin;
267   SetStatus(kStopped);
268 }
269
270
271 //_____________________________________________________________________________
272 Bool_t AliMonitorProcess::CheckForNewFile()
273 {
274 // check whether a new file was registered in alien
275
276   TGridResult* result = fGrid->Ls();
277   Long_t maxDate = -1;
278   Long_t maxTime = -1;
279   TString fileName;
280
281   while (const char* entry = result->Next()) {
282     // entry = host_date_time.root
283     TString entryCopy(entry);
284     char* p = const_cast<char*>(entryCopy.Data());
285     if (!strtok(p, "_") || !p) continue;  // host name
286     char* dateStr = strtok(NULL, "_");
287     if (!dateStr || !p) continue;
288     char* timeStr = strtok(NULL, ".");
289     if (!timeStr || !p) continue;
290     Long_t date = atoi(dateStr);
291     Long_t time = atoi(timeStr);
292
293     if ((date > maxDate) || ((date == maxDate) && (time > maxTime))) {
294       maxDate = date;
295       maxTime = time;
296       fileName = entry;
297     }
298   }
299
300   if (maxDate < 0) return kFALSE;  // no files found
301   if (fLogicalFileName.CompareTo(fileName) == 0) return kFALSE;  // no new file
302
303   fLogicalFileName = fileName;
304   TGridResult* result2 = fGrid->GetPhysicalFileNames(fLogicalFileName.Data());
305   fFileName = result2->Next();
306
307   return kTRUE;
308 }
309
310 //_____________________________________________________________________________
311 Bool_t AliMonitorProcess::ProcessFile()
312 {
313 // loop over all events in the raw data file, run the reconstruction
314 // and fill the monitor histograms
315
316   Int_t nEvents = GetNumberOfEvents(fFileName);
317   if (nEvents <= 0) return kFALSE;
318   Info("ProcessFile", "found %d event(s) in file %s", 
319        nEvents, fFileName.Data());
320 #ifdef ALI_HLT
321   CreateHLT(fFileName);
322   CreateHLTHough(fFileName);
323 #endif
324
325   // loop over the events
326   for (Int_t iEvent = 0; iEvent < nEvents; iEvent++) {
327     SetStatus(kReading);
328     fRunLoader->SetEventNumber(0);
329     AliRawReaderRoot rawReader(fFileName, iEvent);
330     if (fStopping) break;
331     if (rawReader.GetRunNumber() != fRunNumber) {
332       WriteHistos();
333       StartNewRun();
334       fRunNumber = rawReader.GetRunNumber();
335       fEventNumber[0] = rawReader.GetEventId()[0];
336       fEventNumber[1] = rawReader.GetEventId()[1];
337       fSubRunNumber = 0;
338       if (fStopping) break;
339     }
340
341     if (!ReconstructTPC(&rawReader)) return kFALSE;
342     if (fStopping) break;
343     if (!ReconstructITS(&rawReader)) return kFALSE;
344     if (fStopping) break;
345     if (!ReconstructV0s()) return kFALSE;
346     if (fStopping) break;
347     if (!ReconstructHLT(iEvent)) return kFALSE;
348     if (fStopping) break;
349     if (!ReconstructHLTHough(iEvent)) return kFALSE;
350     if (fStopping) break;
351
352     if (fDisplaySocket) fDisplaySocket->Send("new event");
353
354     Info("ProcessFile", "filling histograms...");
355     SetStatus(kFilling);
356     for (Int_t iMonitor = 0; iMonitor < fMonitors.GetEntriesFast(); iMonitor++) {
357       ((AliMonitor*) fMonitors[iMonitor])->FillHistos(fRunLoader, &rawReader);
358       if (fStopping) break;
359     }
360     if (fStopping) break;
361
362     Info("ProcessFile", "updating histograms...");
363     SetStatus(kUpdating);
364     TIterator* iFolder = fTopFolder->GetListOfFolders()->MakeIterator();
365     while (TFolder* folder = (TFolder*) iFolder->Next()) {
366       TIterator* iHisto = folder->GetListOfFolders()->MakeIterator();
367       while (AliMonitorPlot* histo = (AliMonitorPlot*) iHisto->Next()) {
368         histo->Update();
369       }
370       delete iHisto;
371     }
372     delete iFolder;
373     if (fStopping) break;
374
375     Info("ProcessFile", "filling the tree...");
376     fTree->Fill();
377
378     Info("ProcessFile", "broadcasting histograms...");
379     CheckForConnections();
380     BroadcastHistos();
381
382     fNEvents++;
383     if (fStopping) break;
384   }
385
386 #ifdef ALI_HLT
387   delete fHLT;
388   fHLT = NULL;
389 #endif
390
391   return kTRUE;
392 }
393
394 //_____________________________________________________________________________
395 void AliMonitorProcess::Reset()
396 {
397 // write the current histograms to a file and reset them
398
399   if (fSubRunNumber == 0) fSubRunNumber++;
400   WriteHistos();
401   StartNewRun();
402   fSubRunNumber++;
403 }
404
405
406 //_____________________________________________________________________________
407 UInt_t AliMonitorProcess::GetEventPeriodNumber() const
408 {
409 // get the period number from the event id
410
411   return (fEventNumber[1] >> 4);
412 }
413
414 //_____________________________________________________________________________
415 UInt_t AliMonitorProcess::GetEventOrbitNumber() const
416 {
417 // get the orbit number from the event id
418
419   return ((fEventNumber[1] & 0x000F) << 20) + (fEventNumber[0] >> 12);
420 }
421
422 //_____________________________________________________________________________
423 UInt_t AliMonitorProcess::GetEventBunchNumber() const
424 {
425 // get the bunch number from the event id
426
427   return (fEventNumber[0] % 0x0FFF);
428 }
429
430 //_____________________________________________________________________________
431 Int_t AliMonitorProcess::GetNumberOfEvents(const char* fileName) const
432 {
433 // determine the number of events in the given raw data file
434
435   Int_t nEvents = -1;
436
437   TFile* file = TFile::Open(fileName);
438   if (!file || !file->IsOpen()) {
439     Error("GetNumberOfEvents", "could not open file %s", fileName);
440     if (file) delete file;
441     return -1;
442   }
443
444   TTree* tree = (TTree*) file->Get("RAW");
445   if (!tree) {
446     Error("GetNumberOfEvents", "could not find tree with raw data");
447   } else {
448     nEvents = (Int_t) tree->GetEntries();
449   }
450   file->Close();
451   delete file;
452
453   return nEvents;
454 }
455
456 //_____________________________________________________________________________
457 Bool_t AliMonitorProcess::ReconstructTPC(AliRawReader* rawReader)
458 {
459 // find TPC clusters and tracks
460
461   SetStatus(kRecTPC);
462
463   AliLoader* tpcLoader = fRunLoader->GetLoader("TPCLoader");
464   if (!tpcLoader) {
465     Error("ReconstructTPC", "no TPC loader found");
466     return kFALSE;
467   }
468   gSystem->Unlink("TPC.RecPoints.root");
469   gSystem->Unlink("TPC.Tracks.root");
470
471   // cluster finder
472   Info("ReconstructTPC", "reconstructing clusters...");
473   tpcLoader->LoadRecPoints("recreate");
474   AliTPCclustererMI clusterer(fTPCParam);
475   tpcLoader->MakeRecPointsContainer();
476   clusterer.SetOutput(tpcLoader->TreeR());
477   clusterer.Digits2Clusters(rawReader);
478   tpcLoader->WriteRecPoints("OVERWRITE");
479
480   // track finder
481   Info("ReconstructTPC", "reconstructing tracks...");
482   tpcLoader->LoadTracks("recreate");
483   {
484     AliTPCtrackerMI tracker(fTPCParam);
485     tracker.Clusters2Tracks();
486   }
487
488   tpcLoader->UnloadRecPoints();
489   tpcLoader->UnloadTracks();
490   return kTRUE;
491 }
492
493 //_____________________________________________________________________________
494 Bool_t AliMonitorProcess::ReconstructITS(AliRawReader* rawReader)
495 {
496 // find ITS clusters and tracks
497
498   SetStatus(kRecITS);
499
500   AliLoader* itsLoader = fRunLoader->GetLoader("ITSLoader");
501   if (!itsLoader) {
502     Error("ReconstructITS", "no ITS loader found");
503     return kFALSE;
504   }
505   AliLoader* tpcLoader = fRunLoader->GetLoader("TPCLoader");
506   if (!tpcLoader) {
507     Error("ReconstructITS", "no TPC loader found");
508     return kFALSE;
509   }
510   gSystem->Unlink("ITS.RecPoints.root");
511   gSystem->Unlink("ITS.Tracks.root");
512
513   // cluster finder
514   Info("ReconstructITS", "reconstructing clusters...");
515   itsLoader->LoadRecPoints("recreate");
516   AliITSclustererV2 clusterer(fITSgeom);
517   itsLoader->MakeRecPointsContainer();
518   clusterer.Digits2Clusters(rawReader);
519
520   // track finder
521   Info("ReconstructITS", "reconstructing tracks...");
522   itsLoader->LoadTracks("recreate");
523   itsLoader->MakeTracksContainer();
524   tpcLoader->LoadTracks();
525   AliITStrackerV2 tracker(fITSgeom);
526   tracker.LoadClusters(itsLoader->TreeR());
527   tracker.Clusters2Tracks(tpcLoader->TreeT(), itsLoader->TreeT());
528   tracker.UnloadClusters();
529   itsLoader->WriteTracks("OVERWRITE");
530
531   itsLoader->UnloadRecPoints();
532   itsLoader->UnloadTracks();
533   tpcLoader->UnloadTracks();
534   return kTRUE;
535 }
536
537 //_____________________________________________________________________________
538 Bool_t AliMonitorProcess::ReconstructV0s()
539 {
540 // find V0s
541
542   SetStatus(kRecV0s);
543
544   AliITSLoader* itsLoader = (AliITSLoader*) fRunLoader->GetLoader("ITSLoader");
545   if (!itsLoader) {
546     Error("ReconstructV0", "no ITS loader found");
547     return kFALSE;
548   }
549   gSystem->Unlink("ITS.V0s.root");
550
551   // V0 finder
552   Info("ReconstructV0s", "reconstructing V0s...");
553   itsLoader->LoadTracks("read");
554   itsLoader->LoadV0s("recreate");
555   AliV0vertexer vertexer;
556   TTree* tracks = itsLoader->TreeT();
557   if (!tracks) {
558     Error("ReconstructV0s", "no ITS tracks tree found");
559     return kFALSE;
560   }
561   if (!itsLoader->TreeV0()) itsLoader->MakeTree("V0");
562   TTree* v0s = itsLoader->TreeV0();
563   vertexer.Tracks2V0vertices(tracks, v0s);
564   itsLoader->WriteV0s("OVERWRITE");
565
566   itsLoader->UnloadTracks();
567   itsLoader->UnloadV0s();
568   return kTRUE;
569 }
570
571 //_____________________________________________________________________________
572 #ifdef ALI_HLT
573 void AliMonitorProcess::CreateHLT(const char* fileName)
574 {
575
576 // create the HLT (Level3) object
577
578   if (fHLT) delete fHLT;
579
580   char name[256];
581   strcpy(name, fileName);
582   fHLT = new AliLevel3(name);
583   fHLT->Init("./", AliLevel3::kRaw, 1);
584
585   fHLT->SetClusterFinderParam(-1, -1, kTRUE);
586   
587   Int_t phiSegments = 50;
588   Int_t etaSegments = 100;
589   Int_t trackletlength = 3;
590   Int_t tracklength = 20;//40 or 5
591   Int_t rowscopetracklet = 2;
592   Int_t rowscopetrack = 10;
593   Double_t minPtFit = 0;
594   Double_t maxangle = 0.1745;
595   Double_t goodDist = 5;
596   Double_t maxphi = 0.1;
597   Double_t maxeta = 0.1;
598   Double_t hitChi2Cut = 15;//100 or 15
599   Double_t goodHitChi2 = 5;//20 or 5
600   Double_t trackChi2Cut = 10;//50 or 10
601   fHLT->SetTrackerParam(phiSegments, etaSegments, 
602                         trackletlength, tracklength,
603                         rowscopetracklet, rowscopetrack,
604                         minPtFit, maxangle, goodDist, hitChi2Cut,
605                         goodHitChi2, trackChi2Cut, 50, maxphi, maxeta, kTRUE);
606   
607   fHLT->WriteFiles("./hlt/");  
608 }
609
610 //_____________________________________________________________________________
611 void AliMonitorProcess::CreateHLTHough(const char* fileName)
612 {
613
614 // create the HLT Hough transform (L3Hough) object
615
616   if (fHLTHough) delete fHLTHough;
617
618   char name[256];
619   strcpy(name, fileName);
620
621   fHLTHough = new AliL3Hough();
622   fHLTHough->SetTransformerParams(64,64,0.1,30);
623 //  fHLTHough->Init("./", kFALSE, 20, kFALSE,0,name);
624   fHLTHough->SetAddHistograms();
625   fHLTHough->GetMaxFinder()->SetThreshold(55000); // or 14000 ?
626
627 }
628 #endif
629
630 //_____________________________________________________________________________
631 Bool_t AliMonitorProcess::ReconstructHLT(
632 #ifdef ALI_HLT
633   Int_t iEvent
634 #else
635   Int_t /* iEvent */
636 #endif
637 )
638 {
639 // run the HLT cluster and track finder
640
641   SetStatus(kRecHLT);
642
643 #ifndef ALI_HLT
644   Warning("ReconstructHLT", "the code was compiled without HLT support");
645   return kTRUE;
646
647 #else
648   gSystem->Exec("rm -rf hlt");
649   gSystem->MakeDirectory("hlt");
650   if (!fHLT) return kFALSE;
651
652   fHLT->ProcessEvent(0, 35, iEvent);
653
654   // remove the event number from the file names
655   char command[256];
656   sprintf(command, "rename points_%d points hlt/*.raw", iEvent);
657   gSystem->Exec(command);
658   sprintf(command, "rename tracks_tr_%d tracks_tr hlt/*.raw", iEvent);
659   gSystem->Exec(command);
660   sprintf(command, "rename tracks_gl_%d tracks_gl hlt/*.raw", iEvent);
661   gSystem->Exec(command);
662   sprintf(command, "rename tracks_%d tracks hlt/*.raw", iEvent);
663   gSystem->Exec(command);
664   return kTRUE;
665 #endif
666 }
667
668 //_____________________________________________________________________________
669 Bool_t AliMonitorProcess::ReconstructHLTHough(
670 #ifdef ALI_HLT
671   Int_t iEvent
672 #else
673   Int_t /* iEvent */
674 #endif
675 )
676 {
677 // run the HLT Hough transformer
678
679   SetStatus(kRecHLT);
680
681 #ifndef ALI_HLT
682   Warning("ReconstructHLTHough", "the code was compiled without HLT support");
683   return kTRUE;
684
685 #else
686   if (!fHLTHough) return kFALSE;
687
688   //  fHLTHough->Process(0, 35);
689   // Loop over TPC sectors and process the data
690   for(Int_t i=0; i<=35; i++)
691     {
692       fHLTHough->ReadData(i,iEvent);
693       fHLTHough->Transform();
694       //      if(fHLTHough->fAddHistograms)
695       fHLTHough->AddAllHistograms();
696       fHLTHough->FindTrackCandidates();
697       fHLTHough->WriteTracks(i,"./hlt");
698     }
699   return kTRUE;
700 #endif
701 }
702
703 //_____________________________________________________________________________
704 Bool_t AliMonitorProcess::WriteHistos()
705 {
706 // write the monitor tree and the monitor histograms to the file 
707 // "monitor_<run number>[_<sub_run_number>].root"
708 // if at least fNEventsMin events were monitored
709
710   SetStatus(kWriting);
711
712   // rename tree file and create a new one
713   fFile->cd();
714   fTree->Write();
715   fFile->Close();
716   delete fFile;
717
718   char fileName[256];
719   sprintf(fileName, "monitor_tree_%d.root", fRunNumber);
720   if (fSubRunNumber > 0) {
721     sprintf(fileName, "monitor_tree_%d_%d.root", fRunNumber, fSubRunNumber);
722   }
723   if (fNEvents < fNEventsMin) {
724     gSystem->Unlink("monitor_tree.root");
725   } else {
726     gSystem->Rename("monitor_tree.root", fileName);
727   }
728
729   fFile = TFile::Open("monitor_tree.root", "RECREATE");
730   if (!fFile || !fFile->IsOpen()) {
731     Fatal("WriteHistos", "could not open file for tree");
732   }
733   fTree = new TTree("MonitorTree", "tree for monitoring");
734   for (Int_t iMonitor = 0; iMonitor < fMonitors.GetEntriesFast(); iMonitor++) {
735     ((AliMonitor*) fMonitors[iMonitor])->CreateBranches(fTree);
736   }
737   gROOT->cd();
738
739   // write the histograms
740   if (fNEvents < fNEventsMin) return kTRUE;
741
742   if (!fWriteHistoList) {
743     TIterator* iFolder = fTopFolder->GetListOfFolders()->MakeIterator();
744     while (TFolder* folder = (TFolder*) iFolder->Next()) {
745       TIterator* iHisto = folder->GetListOfFolders()->MakeIterator();
746       while (AliMonitorPlot* histo = (AliMonitorPlot*) iHisto->Next()) {
747         histo->ResetList();
748       }
749       delete iHisto;
750     }
751     delete iFolder;
752   }
753
754   Bool_t result = kTRUE;
755   sprintf(fileName, "monitor_%d.root", fRunNumber);
756   if (fSubRunNumber > 0) {
757     sprintf(fileName, "monitor_%d_%d.root", fRunNumber, fSubRunNumber);
758   }
759   TFile* file = TFile::Open(fileName, "recreate");
760   if (!file || !file->IsOpen()) {
761     Error("WriteHistos", "could not open file %s", fileName);
762     result = kFALSE;
763   } else {
764     fTopFolder->Write();
765     file->Close();
766   }
767   if (file) delete file;
768
769   return result;
770 }
771
772 //_____________________________________________________________________________
773 void AliMonitorProcess::StartNewRun()
774 {
775 // reset the histograms for a new run
776
777   SetStatus(kResetting);
778   TIterator* iFolder = fTopFolder->GetListOfFolders()->MakeIterator();
779   while (TFolder* folder = (TFolder*) iFolder->Next()) {
780     TIterator* iHisto = folder->GetListOfFolders()->MakeIterator();
781     while (AliMonitorPlot* histo = (AliMonitorPlot*) iHisto->Next()) {
782       histo->Reset();
783     }
784     delete iHisto;
785   }
786   delete iFolder;
787
788   fNEvents = 0;
789 }
790
791
792 //_____________________________________________________________________________
793 void AliMonitorProcess::CheckForConnections()
794 {
795 // check if new clients want to connect and add them to the list of sockets
796
797   TMessage message(kMESS_OBJECT);
798   message.WriteObject(fTopFolder); 
799   SetStatus(kConnecting);
800
801   TSocket* socket;
802   while ((socket = fServerSocket->Accept()) != (TSocket*)-1) {
803     socket->SetOption(kNoBlock, 1);
804     char socketType[256];
805     if (!socket->Recv(socketType, 255)) continue;
806     if (strcmp(socketType, "client") == 0) {
807       if ((fNEvents == 0) || (socket->Send(message) >= 0)) {
808         fSockets.Add(socket);
809         TInetAddress adr = socket->GetInetAddress();
810         Info("CheckForConnections", "new client:\n %s (%s), port %d\n",
811              adr.GetHostName(), adr.GetHostAddress(), adr.GetPort());
812       }
813     } else if (strcmp(socketType, "display") == 0) {
814       if (fDisplaySocket) {
815         fDisplaySocket->Close();
816         delete fDisplaySocket;
817       }
818       fDisplaySocket = socket;
819       fDisplaySocket->SetOption(kNoBlock, 1);
820       TInetAddress adr = socket->GetInetAddress();
821       Info("CheckForConnections", "new display:\n %s (%s), port %d\n",
822            adr.GetHostName(), adr.GetHostAddress(), adr.GetPort());
823     }
824   }
825
826   for (Int_t iSocket = 0; iSocket < fSockets.GetEntriesFast(); iSocket++) {
827     socket = (TSocket*) fSockets[iSocket];
828     if (!socket) continue;
829     // remove finished client
830     char str[256];
831     if (socket->Recv(str, 255)) {
832       TString socketMessage(str);
833       if(socketMessage.CompareTo("Finished") == 0) {
834         TInetAddress adr = socket->GetInetAddress();
835         Info("CheckForConnections",
836              "disconnect finished client:\n %s (%s), port %d\n",
837              adr.GetHostName(), adr.GetHostAddress(), adr.GetPort());
838         delete fSockets.RemoveAt(iSocket);
839         continue;
840       }
841     }
842     if (!socket->IsValid()) {
843       // remove invalid sockets from the list
844       TInetAddress adr = socket->GetInetAddress();
845       Info("BroadcastHistos", "disconnect client:\n %s (%s), port %d\n",
846            adr.GetHostName(), adr.GetHostAddress(), adr.GetPort());
847       delete fSockets.RemoveAt(iSocket);
848     }
849   }
850   fSockets.Compress();
851 }
852
853 //_____________________________________________________________________________
854 void AliMonitorProcess::BroadcastHistos()
855 {
856 // send the monitor histograms to the clients
857
858   SetStatus(kBroadcasting);
859   TMessage message(kMESS_OBJECT);
860   message.WriteObject(fTopFolder); 
861
862   for (Int_t iSocket = 0; iSocket < fSockets.GetEntriesFast(); iSocket++) {
863     TSocket* socket = (TSocket*) fSockets[iSocket];
864     if (!socket) continue;
865     socket->SetOption(kNoBlock, 0);
866     if (!socket->IsValid() || (socket->Send(message) < 0)) {
867       // remove the socket from the list if there was an error
868       TInetAddress adr = socket->GetInetAddress();
869       Info("BroadcastHistos", "disconnect client:\n %s (%s), port %d\n",
870            adr.GetHostName(), adr.GetHostAddress(), adr.GetPort());
871       delete fSockets.RemoveAt(iSocket);
872     } else {
873       socket->SetOption(kNoBlock, 1);
874     }
875   }
876   fSockets.Compress();
877 }