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