]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MONITOR/AliMonitorProcess.cxx
AliAODStdParticle -> AliAODParticle; AliAODParticle -> AliVAODParticle
[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 #include <TFile.h>
29 #include <TGrid.h>
30 #include <TGridResult.h>
31 #include <TMessage.h>
32 #include <TROOT.h>
33 #include <TServerSocket.h>
34 #include <TSocket.h>
35
36 #include "AliITS.h"
37 #include "AliITSLoader.h"
38 #include "AliITSclustererV2.h"
39 #include "AliITStrackerV2.h"
40 #include "AliLoader.h"
41 #include "AliMonitorHLT.h"
42 #include "AliMonitorHLTHough.h"
43 #include "AliMonitorITS.h"
44 #include "AliMonitorProcess.h"
45 #include "AliMonitorTPC.h"
46 #include "AliMonitorV0s.h"
47 #include "AliRawReaderRoot.h"
48 #include "AliRun.h"
49 #include "AliTPC.h"
50 #include "AliTPCclustererMI.h"
51 #include "AliTPCtrackerMI.h"
52 #include "AliV0vertexer.h"
53
54 #ifdef ALI_HLT
55 #include <AliL3ClusterFitter.h>
56 #include <AliL3DDLDataFileHandler.h>
57 #include <AliL3Fitter.h>
58 #include <AliL3Hough.h>
59 #include <AliL3HoughBaseTransformer.h>
60 #include <AliL3HoughMaxFinder.h>
61 #include <AliL3StandardIncludes.h>
62 #include <AliL3Track.h>
63 #include <AliL3TrackArray.h>
64 #include <AliL3Transform.h>
65 #include <AliL3Vertex.h>
66 #include <AliLevel3.h>
67 #endif
68
69 ClassImp(AliMonitorProcess) 
70
71
72 const Int_t AliMonitorProcess::fgkPort = 9327;
73
74
75 //_____________________________________________________________________________
76 AliMonitorProcess::AliMonitorProcess(
77 #if ROOT_VERSION_CODE <= 199169   // 3.10/01
78                                      const char* /*alienHost*/,
79 #else
80                                      const char* alienHost,
81 #endif
82                                      const char* alienDir,
83                                      const char* fileNameGalice)
84 {
85 // initialize the monitoring process and the monitor histograms
86
87 #if ROOT_VERSION_CODE <= 199169   // 3.10/01
88   fGrid = TGrid::Connect("alien", gSystem->Getenv("USER"));
89 #else
90   fGrid = TGrid::Connect(alienHost, gSystem->Getenv("USER"));
91 #endif
92   if (!fGrid || fGrid->IsZombie() || !fGrid->IsConnected()) {
93     delete fGrid;
94     Fatal("AliMonitorProcess", "could not connect to alien");
95   }
96 #if ROOT_VERSION_CODE <= 199169   // 3.10/01
97   fGrid->cd(alienDir);
98 #else
99   fAlienDir = alienDir;
100 #endif
101   fLogicalFileName = "";
102   fFileName = "";
103
104   fRunLoader = AliRunLoader::Open(fileNameGalice);
105   if (!fRunLoader) Fatal("AliMonitorProcess", 
106                          "could not get run loader from file %s", 
107                          fileNameGalice);
108
109   fRunLoader->CdGAFile();
110   fTPCParam = AliTPC::LoadTPCParam(gFile);
111   if (!fTPCParam) Fatal("AliMonitorProcess", "could not load TPC parameters");
112
113   fRunLoader->LoadgAlice();
114   gAlice = fRunLoader->GetAliRun();
115   if (!gAlice) Fatal("AliMonitorProcess", "no gAlice object found");
116   AliITS* its = (AliITS*) gAlice->GetModule("ITS");
117   if (!its) Fatal("AliMonitorProcess", "no ITS detector found");
118   fITSgeom = its->GetITSgeom();
119   if (!fITSgeom) Fatal("AliMonitorProcess", "could not load ITS geometry");
120
121 #ifdef ALI_HLT
122 // Init TPC parameters for HLT
123   Bool_t isinit=AliL3Transform::Init(const_cast<char*>(fileNameGalice),kTRUE);
124   if(!isinit){
125     cerr << "Could not create transform settings, please check log for error messages!" << endl;
126     return;
127   }
128 #endif
129
130   fRunNumber = 0;
131   fSubRunNumber = 0;
132   fNEvents = 0;
133   fNEventsMin = 1;
134   fWriteHistoList = kFALSE;
135
136   fTopFolder = new TFolder("Monitor", "monitor histograms");
137   fTopFolder->SetOwner(kTRUE);
138
139   fMonitors.Add(new AliMonitorTPC(fTPCParam));
140   fMonitors.Add(new AliMonitorITS(fITSgeom));
141   fMonitors.Add(new AliMonitorV0s);
142 #ifdef ALI_HLT
143   fMonitors.Add(new AliMonitorHLT(fTPCParam));
144   fMonitors.Add(new AliMonitorHLTHough(fTPCParam));
145 #endif
146
147   for (Int_t iMonitor = 0; iMonitor < fMonitors.GetEntriesFast(); iMonitor++) {
148     ((AliMonitor*) fMonitors[iMonitor])->CreateHistos(fTopFolder);
149   }
150
151   TIterator* iFolder = fTopFolder->GetListOfFolders()->MakeIterator();
152   while (TFolder* folder = (TFolder*) iFolder->Next()) folder->SetOwner(kTRUE);
153   delete iFolder;
154
155   fFile = TFile::Open("monitor_tree.root", "RECREATE");
156   if (!fFile || !fFile->IsOpen()) {
157     Fatal("AliMonitorProcess", "could not open file for tree");
158   }
159   fTree = new TTree("MonitorTree", "tree for monitoring");
160   for (Int_t iMonitor = 0; iMonitor < fMonitors.GetEntriesFast(); iMonitor++) {
161     ((AliMonitor*) fMonitors[iMonitor])->CreateBranches(fTree);
162   }
163   gROOT->cd();
164
165   fServerSocket = new TServerSocket(fgkPort, kTRUE);
166   fServerSocket->SetOption(kNoBlock, 1);
167   fDisplaySocket = NULL;
168   CheckForConnections();
169 #ifdef ALI_HLT
170   fHLT = NULL;
171 #endif
172
173   SetStatus(kStopped);
174   fStopping = kFALSE;
175
176   fInterruptHandler = new AliMonitorInterruptHandler(this);
177   gSystem->AddSignalHandler(fInterruptHandler);
178 }
179
180 //_____________________________________________________________________________
181 AliMonitorProcess::AliMonitorProcess(const AliMonitorProcess& process) :
182   TObject(process)
183 {
184   Fatal("AliMonitorProcess", "copy constructor not implemented");
185 }
186
187 //_____________________________________________________________________________
188 AliMonitorProcess& AliMonitorProcess::operator = (const AliMonitorProcess& 
189                                                   /*process*/)
190 {
191   Fatal("operator =", "assignment operator not implemented");
192   return *this;
193 }
194
195 //_____________________________________________________________________________
196 AliMonitorProcess::~AliMonitorProcess()
197 {
198 // clean up
199
200   fMonitors.Delete();
201   delete fTopFolder;
202   delete fRunLoader;
203
204   delete fServerSocket;
205   fSockets.Delete();
206   delete fDisplaySocket;
207
208 #if ROOT_VERSION_CODE <= 199169   // 3.10/01
209   fGrid->Close();
210 #endif
211   delete fGrid;
212
213   fFile->Close();
214   delete fFile;
215   gSystem->Unlink("monitor_tree.root");
216
217 #ifdef ALI_HLT
218   delete fHLT;
219   delete fHLTHough;
220 #endif
221
222   gSystem->RemoveSignalHandler(fInterruptHandler);
223   delete fInterruptHandler;
224 }
225
226
227 //_____________________________________________________________________________
228 const char* AliMonitorProcess::GetRevision()
229 {
230   return "$Revision$";
231 }
232
233
234 //_____________________________________________________________________________
235 void AliMonitorProcess::SetStatus(EStatus status)
236 {
237 // set the current status and process system events
238
239   fStatus = status;
240   gSystem->ProcessEvents();
241 }
242
243
244 //_____________________________________________________________________________
245 void AliMonitorProcess::Run()
246 {
247 // run the monitor process: 
248 //  check for a raw data file, process the raw data file and delete it
249
250   fStopping = kFALSE;
251
252   while (!fStopping) {
253     SetStatus(kWaiting);
254     while (!CheckForNewFile()) {
255       CheckForConnections();
256       SetStatus(kWaiting);
257       if (fStopping) break;
258       gSystem->Sleep(10);
259     }
260     if (fStopping) break;
261
262     ProcessFile();
263   }
264
265   WriteHistos();
266
267   fStopping = kFALSE;
268   SetStatus(kStopped);
269 }
270
271
272 //_____________________________________________________________________________
273 void AliMonitorProcess::Stop()
274 {
275 // set the fStopping flag to terminate the monitor process after the current
276 // event was processed
277
278   if (GetStatus() != kStopped) fStopping = kTRUE;
279 }
280
281
282 //_____________________________________________________________________________
283 void AliMonitorProcess::ProcessFile(const char* fileName)
284 {
285 // create a file with monitor histograms for a single file
286
287   if (GetStatus() != kStopped) {
288     Error("ProcessFile", "ProcessFile can not be called"
289           " while the monitor process is running");
290     return;
291   }
292
293   fFileName = fileName;
294   Int_t nEventMin = fNEventsMin;
295   fNEventsMin = 1;
296   ProcessFile();
297   WriteHistos();
298   fNEventsMin = nEventMin;
299   SetStatus(kStopped);
300 }
301
302
303 //_____________________________________________________________________________
304 Bool_t AliMonitorProcess::CheckForNewFile()
305 {
306 // check whether a new file was registered in alien
307
308 #if ROOT_VERSION_CODE <= 199169   // 3.10/01
309   TGridResult* result = fGrid->Ls();
310 #else
311   TDatime datime;
312   char dirName[256];
313   sprintf(dirName, "%s/adc-%d", fAlienDir.Data(), datime.GetDate());
314   char findName[256];
315   sprintf(findName, "*.root");
316   Grid_ResultHandle_t handle = fGrid->Find(dirName, findName);
317   if (!handle) {
318     Error("CheckForNewFile", "could not open alien directory %s", 
319           dirName);
320     return kFALSE;
321   }
322   TGridResult* result = fGrid->CreateGridResult(handle);
323 #endif
324   Long_t maxDate = -1;
325   Long_t maxTime = -1;
326   TString fileName;
327
328 #if ROOT_VERSION_CODE <= 199169   // 3.10/01
329   while (const char* entry = result->Next()) {
330 #else
331   while (Grid_Result_t* resultEntry = result->Next()) {
332     const char* entry = resultEntry->name.c_str();
333 #endif
334     if (strrchr(entry, '/')) entry = strrchr(entry, '/')+1;
335     // entry = host_date_time.root
336     TString entryCopy(entry);
337     char* p = const_cast<char*>(entryCopy.Data());
338     if (!strtok(p, "_") || !p) continue;  // host name
339     char* dateStr = strtok(NULL, "_");
340     if (!dateStr || !p) continue;
341     char* timeStr = strtok(NULL, ".");
342     if (!timeStr || !p) continue;
343     Long_t date = atoi(dateStr);
344     Long_t time = atoi(timeStr);
345
346     if ((date > maxDate) || ((date == maxDate) && (time > maxTime))) {
347       maxDate = date;
348       maxTime = time;
349       fileName = entry;
350     }
351   }
352
353   delete result;
354   if (maxDate < 0) return kFALSE;  // no files found
355   if (fLogicalFileName.CompareTo(fileName) == 0) return kFALSE;  // no new file
356
357   fLogicalFileName = fileName;
358 #if ROOT_VERSION_CODE <= 199169   // 3.10/01
359   result = fGrid->GetPhysicalFileNames(fLogicalFileName.Data());
360   fFileName = result->Next();
361 #else
362   fileName = dirName + ("/" + fLogicalFileName);
363   handle = fGrid->GetPhysicalFileNames(fileName.Data());
364   if (!handle) {
365     Error("CheckForNewFile", "could not get physical file names for %s", 
366           fileName.Data());
367     return kFALSE;
368   }
369   result = fGrid->CreateGridResult(handle);
370   result->Reset();
371   Grid_Result_t* resultEntry = result->Next();
372   if (!resultEntry) {
373     Error("CheckForNewFile", "could not get physical file names for %s", 
374           fileName.Data());
375     return kFALSE;
376   }
377   fFileName = resultEntry->name2.c_str();
378   fFileName.ReplaceAll("castor:/", "rfio:/");
379 #endif
380   delete result;
381
382   return kTRUE;
383 }
384
385 //_____________________________________________________________________________
386 Bool_t AliMonitorProcess::ProcessFile()
387 {
388 // loop over all events in the raw data file, run the reconstruction
389 // and fill the monitor histograms
390
391   Int_t nEvents = GetNumberOfEvents(fFileName);
392   if (nEvents <= 0) return kFALSE;
393   Info("ProcessFile", "found %d event(s) in file %s", 
394        nEvents, fFileName.Data());
395 #ifdef ALI_HLT
396   CreateHLT(fFileName);
397   CreateHLTHough(fFileName);
398 #endif
399
400   // loop over the events
401   for (Int_t iEvent = 0; iEvent < nEvents; iEvent++) {
402     CheckForConnections();
403     SetStatus(kReading);
404     fRunLoader->SetEventNumber(0);
405     AliRawReaderRoot rawReader(fFileName, iEvent);
406     if (fStopping) break;
407     if (rawReader.GetRunNumber() != fRunNumber) {
408       WriteHistos();
409       StartNewRun();
410       fRunNumber = rawReader.GetRunNumber();
411       fEventNumber[0] = rawReader.GetEventId()[0];
412       fEventNumber[1] = rawReader.GetEventId()[1];
413       fSubRunNumber = 0;
414       if (fStopping) break;
415     }
416
417     // monitor only central physics events
418     if (rawReader.GetType() != 7) continue;
419     if ((rawReader.GetAttributes()[0] & 0x02) == 0) continue;
420     Info("ProcessFile", "run: %d  event: %d %d\n", rawReader.GetRunNumber(), 
421          rawReader.GetEventId()[0], rawReader.GetEventId()[1]);
422
423     CheckForConnections();
424     if (!ReconstructTPC(&rawReader)) return kFALSE;
425     if (fStopping) break;
426     CheckForConnections();
427     if (!ReconstructITS(&rawReader)) return kFALSE;
428     if (fStopping) break;
429     CheckForConnections();
430     if (!ReconstructV0s()) return kFALSE;
431     if (fStopping) break;
432     CheckForConnections();
433     if (!ReconstructHLT(iEvent)) return kFALSE;
434     if (fStopping) break;
435     CheckForConnections();
436     if (!ReconstructHLTHough(iEvent)) return kFALSE;
437     if (fStopping) break;
438
439     if (fDisplaySocket) fDisplaySocket->Send("new event");
440
441     Info("ProcessFile", "filling histograms...");
442     for (Int_t iMonitor = 0; iMonitor < fMonitors.GetEntriesFast(); iMonitor++) {
443       CheckForConnections();
444       SetStatus(kFilling);
445       ((AliMonitor*) fMonitors[iMonitor])->FillHistos(fRunLoader, &rawReader);
446       if (fStopping) break;
447     }
448     if (fStopping) break;
449
450     Info("ProcessFile", "updating histograms...");
451     CheckForConnections();
452     SetStatus(kUpdating);
453     TIterator* iFolder = fTopFolder->GetListOfFolders()->MakeIterator();
454     while (TFolder* folder = (TFolder*) iFolder->Next()) {
455       TIterator* iHisto = folder->GetListOfFolders()->MakeIterator();
456       while (AliMonitorPlot* histo = (AliMonitorPlot*) iHisto->Next()) {
457         histo->Update();
458       }
459       delete iHisto;
460     }
461     delete iFolder;
462     if (fStopping) break;
463
464     Info("ProcessFile", "filling the tree...");
465     fTree->Fill();
466
467     Info("ProcessFile", "broadcasting histograms...");
468     CheckForConnections();
469     BroadcastHistos();
470
471     fNEvents++;
472     if (fStopping) break;
473   }
474
475 #ifdef ALI_HLT
476   delete fHLT;
477   fHLT = NULL;
478 #endif
479
480   return kTRUE;
481 }
482
483 //_____________________________________________________________________________
484 void AliMonitorProcess::Reset()
485 {
486 // write the current histograms to a file and reset them
487
488   if (fSubRunNumber == 0) fSubRunNumber++;
489   WriteHistos();
490   StartNewRun();
491   fSubRunNumber++;
492 }
493
494
495 //_____________________________________________________________________________
496 UInt_t AliMonitorProcess::GetEventPeriodNumber() const
497 {
498 // get the period number from the event id
499
500   return (fEventNumber[1] >> 4);
501 }
502
503 //_____________________________________________________________________________
504 UInt_t AliMonitorProcess::GetEventOrbitNumber() const
505 {
506 // get the orbit number from the event id
507
508   return ((fEventNumber[1] & 0x000F) << 20) + (fEventNumber[0] >> 12);
509 }
510
511 //_____________________________________________________________________________
512 UInt_t AliMonitorProcess::GetEventBunchNumber() const
513 {
514 // get the bunch number from the event id
515
516   return (fEventNumber[0] % 0x0FFF);
517 }
518
519 //_____________________________________________________________________________
520 Int_t AliMonitorProcess::GetNumberOfEvents(const char* fileName) const
521 {
522 // determine the number of events in the given raw data file
523
524   Int_t nEvents = -1;
525
526   TFile* file = TFile::Open(fileName);
527   if (!file || !file->IsOpen()) {
528     Error("GetNumberOfEvents", "could not open file %s", fileName);
529     if (file) delete file;
530     return -1;
531   }
532
533   TTree* tree = (TTree*) file->Get("RAW");
534   if (!tree) {
535     Error("GetNumberOfEvents", "could not find tree with raw data");
536   } else {
537     nEvents = (Int_t) tree->GetEntries();
538   }
539   file->Close();
540   delete file;
541
542   return nEvents;
543 }
544
545 //_____________________________________________________________________________
546 Bool_t AliMonitorProcess::ReconstructTPC(AliRawReader* rawReader)
547 {
548 // find TPC clusters and tracks
549
550   SetStatus(kRecTPC);
551
552   AliLoader* tpcLoader = fRunLoader->GetLoader("TPCLoader");
553   if (!tpcLoader) {
554     Error("ReconstructTPC", "no TPC loader found");
555     return kFALSE;
556   }
557   gSystem->Unlink("TPC.RecPoints.root");
558   gSystem->Unlink("TPC.Tracks.root");
559
560   // cluster finder
561   Info("ReconstructTPC", "reconstructing clusters...");
562   tpcLoader->LoadRecPoints("recreate");
563   AliTPCclustererMI clusterer(fTPCParam);
564   tpcLoader->MakeRecPointsContainer();
565   clusterer.SetOutput(tpcLoader->TreeR());
566   clusterer.Digits2Clusters(rawReader);
567   tpcLoader->WriteRecPoints("OVERWRITE");
568
569   // track finder
570   Info("ReconstructTPC", "reconstructing tracks...");
571   tpcLoader->LoadTracks("recreate");
572   {
573     AliTPCtrackerMI tracker(fTPCParam);
574     tracker.SetIO();
575     tracker.LoadClusters();
576     tracker.Clusters2Tracks();
577     tracker.WriteTracks();
578     tracker.UnloadClusters();
579     tpcLoader->WriteTracks("OVERWRITE");
580   }
581
582   tpcLoader->UnloadRecPoints();
583   tpcLoader->UnloadTracks();
584   return kTRUE;
585 }
586
587 //_____________________________________________________________________________
588 Bool_t AliMonitorProcess::ReconstructITS(AliRawReader* rawReader)
589 {
590 // find ITS clusters and tracks
591
592   SetStatus(kRecITS);
593
594   AliLoader* itsLoader = fRunLoader->GetLoader("ITSLoader");
595   if (!itsLoader) {
596     Error("ReconstructITS", "no ITS loader found");
597     return kFALSE;
598   }
599   AliLoader* tpcLoader = fRunLoader->GetLoader("TPCLoader");
600   if (!tpcLoader) {
601     Error("ReconstructITS", "no TPC loader found");
602     return kFALSE;
603   }
604   gSystem->Unlink("ITS.RecPoints.root");
605   gSystem->Unlink("ITS.Tracks.root");
606
607   // cluster finder
608   Info("ReconstructITS", "reconstructing clusters...");
609   itsLoader->LoadRecPoints("recreate");
610   AliITSclustererV2 clusterer(fITSgeom);
611   itsLoader->MakeRecPointsContainer();
612   clusterer.Digits2Clusters(rawReader);
613
614   // track finder
615   Info("ReconstructITS", "reconstructing tracks...");
616   itsLoader->LoadTracks("recreate");
617   itsLoader->MakeTracksContainer();
618   tpcLoader->LoadTracks();
619   AliITStrackerV2 tracker(fITSgeom);
620   tracker.LoadClusters(itsLoader->TreeR());
621   tracker.Clusters2Tracks(tpcLoader->TreeT(), itsLoader->TreeT());
622   tracker.UnloadClusters();
623   itsLoader->WriteTracks("OVERWRITE");
624
625   itsLoader->UnloadRecPoints();
626   itsLoader->UnloadTracks();
627   tpcLoader->UnloadTracks();
628   return kTRUE;
629 }
630
631 //_____________________________________________________________________________
632 Bool_t AliMonitorProcess::ReconstructV0s()
633 {
634 // find V0s
635
636   SetStatus(kRecV0s);
637
638   AliITSLoader* itsLoader = (AliITSLoader*) fRunLoader->GetLoader("ITSLoader");
639   if (!itsLoader) {
640     Error("ReconstructV0", "no ITS loader found");
641     return kFALSE;
642   }
643   gSystem->Unlink("ITS.V0s.root");
644
645   // V0 finder
646   Info("ReconstructV0s", "reconstructing V0s...");
647   itsLoader->LoadTracks("read");
648   itsLoader->LoadV0s("recreate");
649   AliV0vertexer vertexer;
650   TTree* tracks = itsLoader->TreeT();
651   if (!tracks) {
652     Error("ReconstructV0s", "no ITS tracks tree found");
653     return kFALSE;
654   }
655   if (!itsLoader->TreeV0()) itsLoader->MakeTree("V0");
656   TTree* v0s = itsLoader->TreeV0();
657   vertexer.Tracks2V0vertices(tracks, v0s);
658   itsLoader->WriteV0s("OVERWRITE");
659
660   itsLoader->UnloadTracks();
661   itsLoader->UnloadV0s();
662   return kTRUE;
663 }
664
665 //_____________________________________________________________________________
666 #ifdef ALI_HLT
667 void AliMonitorProcess::CreateHLT(const char* fileName)
668 {
669
670 // create the HLT (Level3) object
671
672   if (fHLT) delete fHLT;
673
674   char name[256];
675   strcpy(name, fileName);
676   fHLT = new AliLevel3(name);
677   fHLT->Init("./", AliLevel3::kRaw, 1);
678
679   fHLT->SetClusterFinderParam(-1, -1, kTRUE);
680   
681   Int_t phiSegments = 50;
682   Int_t etaSegments = 100;
683   Int_t trackletlength = 3;
684   Int_t tracklength = 20;//40 or 5
685   Int_t rowscopetracklet = 2;
686   Int_t rowscopetrack = 10;
687   Double_t minPtFit = 0;
688   Double_t maxangle = 0.1745;
689   Double_t goodDist = 5;
690   Double_t maxphi = 0.1;
691   Double_t maxeta = 0.1;
692   Double_t hitChi2Cut = 15;//100 or 15
693   Double_t goodHitChi2 = 5;//20 or 5
694   Double_t trackChi2Cut = 10;//50 or 10
695   fHLT->SetTrackerParam(phiSegments, etaSegments, 
696                         trackletlength, tracklength,
697                         rowscopetracklet, rowscopetrack,
698                         minPtFit, maxangle, goodDist, hitChi2Cut,
699                         goodHitChi2, trackChi2Cut, 50, maxphi, maxeta, kTRUE);
700   
701   fHLT->WriteFiles("./hlt/");  
702 }
703
704 //_____________________________________________________________________________
705 void AliMonitorProcess::CreateHLTHough(const char* fileName)
706 {
707
708 // create the HLT Hough transform (L3Hough) object
709
710   if (fHLTHough) delete fHLTHough;
711
712   char name[256];
713   strcpy(name, fileName);
714
715   fHLTHough = new AliL3Hough();
716   fHLTHough->SetThreshold(4);
717   fHLTHough->SetTransformerParams(140,150,0.5,-1);
718   fHLTHough->SetPeakThreshold(9000,-1);// or 6000
719   fHLTHough->Init("./", kFALSE, 50, kFALSE,0,name);
720   fHLTHough->SetAddHistograms();
721   //  fHLTHough->GetMaxFinder()->SetThreshold(14000);
722
723 }
724 #endif
725
726 //_____________________________________________________________________________
727 Bool_t AliMonitorProcess::ReconstructHLT(
728 #ifdef ALI_HLT
729   Int_t iEvent
730 #else
731   Int_t /* iEvent */
732 #endif
733 )
734 {
735 // run the HLT cluster and track finder
736
737   SetStatus(kRecHLT);
738
739 #ifndef ALI_HLT
740   Warning("ReconstructHLT", "the code was compiled without HLT support");
741   return kTRUE;
742
743 #else
744   gSystem->Exec("rm -rf hlt");
745   gSystem->MakeDirectory("hlt");
746   if (!fHLT) return kFALSE;
747
748   fHLT->ProcessEvent(0, 35, iEvent);
749
750   // remove the event number from the file names
751   char command[256];
752   sprintf(command, "rename points_%d points hlt/*.raw", iEvent);
753   gSystem->Exec(command);
754   sprintf(command, "rename tracks_tr_%d tracks_tr hlt/*.raw", iEvent);
755   gSystem->Exec(command);
756   sprintf(command, "rename tracks_gl_%d tracks_gl hlt/*.raw", iEvent);
757   gSystem->Exec(command);
758   sprintf(command, "rename tracks_%d tracks hlt/*.raw", iEvent);
759   gSystem->Exec(command);
760   return kTRUE;
761 #endif
762 }
763
764 //_____________________________________________________________________________
765 Bool_t AliMonitorProcess::ReconstructHLTHough(
766 #ifdef ALI_HLT
767   Int_t iEvent
768 #else
769   Int_t /* iEvent */
770 #endif
771 )
772 {
773 // run the HLT Hough transformer
774
775   SetStatus(kRecHLT);
776
777 #ifndef ALI_HLT
778   Warning("ReconstructHLTHough", "the code was compiled without HLT support");
779   return kTRUE;
780
781 #else
782   gSystem->Exec("rm -rf hlt/hough");
783   gSystem->MakeDirectory("hlt/hough");
784   gSystem->Exec("rm -rf hlt/fitter");
785   gSystem->MakeDirectory("hlt/fitter");
786   if (!fHLTHough) return kFALSE;
787
788   //  fHLTHough->Process(0, 35);
789   // Loop over TPC sectors and process the data
790   for(Int_t i=0; i<=35; i++)
791     {
792       fHLTHough->ReadData(i,iEvent);
793       fHLTHough->Transform();
794       //      if(fHLTHough->fAddHistograms)
795       fHLTHough->AddAllHistograms();
796       fHLTHough->FindTrackCandidates();
797       fHLTHough->AddTracks();
798     }
799   fHLTHough->WriteTracks("./hlt/hough");
800
801   // Run cluster fitter
802   AliL3ClusterFitter *fitter = new AliL3ClusterFitter("./hlt");
803
804   // Set debug flag for the cluster fitter
805   //  fitter->Debug();
806
807   // Setting fitter parameters
808   fitter->SetInnerWidthFactor(1,1.5);
809   fitter->SetOuterWidthFactor(1,1.5);
810   fitter->SetNmaxOverlaps(5);
811   
812   //fitter->SetChiSqMax(5,kFALSE); //isolated clusters
813   fitter->SetChiSqMax(5,kTRUE);  //overlapping clusters
814
815   Int_t rowrange[2] = {0,AliL3Transform::GetNRows()-1};
816
817   // Takes input from global hough tracks produced by HT
818   fitter->LoadSeeds(rowrange,kFALSE,iEvent);
819
820   UInt_t ndigits;
821
822   for(Int_t islice = 0; islice <= 35; islice++)
823     {
824       for(Int_t ipatch = 0; ipatch < AliL3Transform::GetNPatches(); ipatch++)
825         {
826           // Read digits
827           fHLTHough->GetMemHandler(ipatch)->Free();
828           fHLTHough->GetMemHandler(ipatch)->Init(islice,ipatch);
829           AliL3DigitRowData *digits = (AliL3DigitRowData *)fHLTHough->GetMemHandler(ipatch)->AliAltroDigits2Memory(ndigits,iEvent);
830
831           fitter->Init(islice,ipatch);
832           fitter->SetInputData(digits);
833           fitter->FindClusters();
834           fitter->WriteClusters();
835         }
836     }
837
838   // Refit of the clusters
839   AliL3Vertex vertex;
840   //The seeds are the input tracks from circle HT
841   AliL3TrackArray *tracks = fitter->GetSeeds();
842   AliL3Fitter *ft = new AliL3Fitter(&vertex,1);
843
844   ft->LoadClusters("./hlt/fitter/",iEvent,kFALSE);
845   for(Int_t i=0; i<tracks->GetNTracks(); i++)
846     {
847       AliL3Track *track = tracks->GetCheckedTrack(i);
848       if(!track) continue;
849       if(track->GetNHits() < 20) continue;
850       ft->SortTrackClusters(track);
851       ft->FitHelix(track);
852       track->UpdateToFirstPoint();
853     }
854   delete ft;
855         
856   //Write the final tracks
857   fitter->WriteTracks(20);
858
859   delete fitter;
860
861   // remove the event number from the file names
862   char command[256];
863   sprintf(command, "rename tracks_%d tracks hlt/hough/*.raw", iEvent);
864   gSystem->Exec(command);
865   sprintf(command, "rename tracks_%d tracks hlt/fitter/*.raw", iEvent);
866   gSystem->Exec(command);
867   sprintf(command, "rename points_%d points hlt/fitter/*.raw", iEvent);
868   gSystem->Exec(command);
869   return kTRUE;
870 #endif
871 }
872
873 //_____________________________________________________________________________
874 Bool_t AliMonitorProcess::WriteHistos()
875 {
876 // write the monitor tree and the monitor histograms to the file 
877 // "monitor_<run number>[_<sub_run_number>].root"
878 // if at least fNEventsMin events were monitored
879
880   SetStatus(kWriting);
881
882   // rename tree file and create a new one
883   fFile->cd();
884   fTree->Write();
885   fFile->Close();
886   delete fFile;
887
888   char fileName[256];
889   sprintf(fileName, "monitor_tree_%d.root", fRunNumber);
890   if (fSubRunNumber > 0) {
891     sprintf(fileName, "monitor_tree_%d_%d.root", fRunNumber, fSubRunNumber);
892   }
893   if (fNEvents < fNEventsMin) {
894     gSystem->Unlink("monitor_tree.root");
895   } else {
896     gSystem->Rename("monitor_tree.root", fileName);
897   }
898
899   fFile = TFile::Open("monitor_tree.root", "RECREATE");
900   if (!fFile || !fFile->IsOpen()) {
901     Fatal("WriteHistos", "could not open file for tree");
902   }
903   fTree = new TTree("MonitorTree", "tree for monitoring");
904   for (Int_t iMonitor = 0; iMonitor < fMonitors.GetEntriesFast(); iMonitor++) {
905     ((AliMonitor*) fMonitors[iMonitor])->CreateBranches(fTree);
906   }
907   gROOT->cd();
908
909   // write the histograms
910   if (fNEvents < fNEventsMin) return kTRUE;
911
912   if (!fWriteHistoList) {
913     TIterator* iFolder = fTopFolder->GetListOfFolders()->MakeIterator();
914     while (TFolder* folder = (TFolder*) iFolder->Next()) {
915       TIterator* iHisto = folder->GetListOfFolders()->MakeIterator();
916       while (AliMonitorPlot* histo = (AliMonitorPlot*) iHisto->Next()) {
917         histo->ResetList();
918       }
919       delete iHisto;
920     }
921     delete iFolder;
922   }
923
924   Bool_t result = kTRUE;
925   sprintf(fileName, "monitor_%d.root", fRunNumber);
926   if (fSubRunNumber > 0) {
927     sprintf(fileName, "monitor_%d_%d.root", fRunNumber, fSubRunNumber);
928   }
929   TFile* file = TFile::Open(fileName, "recreate");
930   if (!file || !file->IsOpen()) {
931     Error("WriteHistos", "could not open file %s", fileName);
932     result = kFALSE;
933   } else {
934     fTopFolder->Write();
935     file->Close();
936   }
937   if (file) delete file;
938
939   return result;
940 }
941
942 //_____________________________________________________________________________
943 void AliMonitorProcess::StartNewRun()
944 {
945 // reset the histograms for a new run
946
947   SetStatus(kResetting);
948   TIterator* iFolder = fTopFolder->GetListOfFolders()->MakeIterator();
949   while (TFolder* folder = (TFolder*) iFolder->Next()) {
950     TIterator* iHisto = folder->GetListOfFolders()->MakeIterator();
951     while (AliMonitorPlot* histo = (AliMonitorPlot*) iHisto->Next()) {
952       histo->Reset();
953     }
954     delete iHisto;
955   }
956   delete iFolder;
957
958   fNEvents = 0;
959 }
960
961
962 //_____________________________________________________________________________
963 void AliMonitorProcess::CheckForConnections()
964 {
965 // check if new clients want to connect and add them to the list of sockets
966
967   TSocket* socket;
968   while ((socket = fServerSocket->Accept()) != (TSocket*)-1) {
969     socket->SetOption(kNoBlock, 1);
970     char socketType[256];
971     if (socket->Recv(socketType, 255) <= 0) {
972       gSystem->Sleep(1000);
973       if (socket->Recv(socketType, 255) <= 0) {
974         TInetAddress adr = socket->GetInetAddress();
975         Error("CheckForConnections", "no socket type received - "
976               "disconnect client:\n %s (%s), port %d\n",
977               adr.GetHostName(), adr.GetHostAddress(), adr.GetPort());
978         delete socket;
979         continue;
980       }
981     }
982     if (strcmp(socketType, "client") == 0) {
983       fSockets.Add(socket);
984       TInetAddress adr = socket->GetInetAddress();
985       Info("CheckForConnections", "new client:\n %s (%s), port %d\n",
986            adr.GetHostName(), adr.GetHostAddress(), adr.GetPort());
987       if (fNEvents > 0) BroadcastHistos(socket);
988     } else if (strcmp(socketType, "display") == 0) {
989       if (fDisplaySocket) {
990         fDisplaySocket->Close();
991         delete fDisplaySocket;
992       }
993       fDisplaySocket = socket;
994       TInetAddress adr = socket->GetInetAddress();
995       Info("CheckForConnections", "new display:\n %s (%s), port %d\n",
996            adr.GetHostName(), adr.GetHostAddress(), adr.GetPort());
997     } else {
998       TInetAddress adr = socket->GetInetAddress();
999       Error("CheckForConnections", "unknown socket type %s - "
1000             "disconnect client:\n %s (%s), port %d\n", socketType,
1001             adr.GetHostName(), adr.GetHostAddress(), adr.GetPort());
1002       delete socket;
1003       continue;
1004     }
1005   }
1006
1007   // remove finished or invalid clients
1008   for (Int_t iSocket = 0; iSocket < fSockets.GetEntriesFast(); iSocket++) {
1009     socket = (TSocket*) fSockets[iSocket];
1010     if (!socket) continue;
1011     char controlMessage[256];
1012     if (socket->Recv(controlMessage, 255)) {
1013       if (strcmp(controlMessage, "disconnect") == 0) {
1014         TInetAddress adr = socket->GetInetAddress();
1015         Info("CheckForConnections",
1016              "disconnect client:\n %s (%s), port %d\n",
1017              adr.GetHostName(), adr.GetHostAddress(), adr.GetPort());
1018         delete fSockets.RemoveAt(iSocket);
1019         continue;
1020       }
1021     }
1022     if (!socket->IsValid()) {
1023       // remove invalid sockets from the list
1024       TInetAddress adr = socket->GetInetAddress();
1025       Error("CheckForConnections", 
1026             "disconnect invalid client:\n %s (%s), port %d\n",
1027             adr.GetHostName(), adr.GetHostAddress(), adr.GetPort());
1028       delete fSockets.RemoveAt(iSocket);
1029     }
1030   }
1031   fSockets.Compress();
1032 }
1033
1034 //_____________________________________________________________________________
1035 void AliMonitorProcess::BroadcastHistos(TSocket* toSocket)
1036 {
1037 // send the monitor histograms to the clients
1038
1039   SetStatus(kBroadcasting);
1040   TMessage message(kMESS_OBJECT);
1041   message.WriteObject(fTopFolder); 
1042
1043   for (Int_t iSocket = 0; iSocket < fSockets.GetEntriesFast(); iSocket++) {
1044     TSocket* socket = (TSocket*) fSockets[iSocket];
1045     if (!socket) continue;
1046     if (toSocket && (socket != toSocket)) continue;
1047
1048     // send control message
1049     if (!socket->IsValid() || (socket->Send("histograms") <= 0)) {
1050       TInetAddress adr = socket->GetInetAddress();
1051       Error("BroadcastHistos", "connection to client failed - "
1052             "disconnect client:\n %s (%s), port %d\n",
1053             adr.GetHostName(), adr.GetHostAddress(), adr.GetPort());
1054       delete fSockets.RemoveAt(iSocket);
1055     }
1056
1057     // receive control message
1058     char controlMessage[256];
1059     Int_t result = socket->Recv(controlMessage, 255);
1060     if (result <= 0) {
1061       gSystem->Sleep(1000);  // wait one second and try again
1062       result = socket->Recv(controlMessage, 255);
1063     }
1064     if (result <= 0) {
1065       TInetAddress adr = socket->GetInetAddress();
1066       Error("BroadcastHistos", "no response from client - "
1067             "disconnect client:\n %s (%s), port %d\n",
1068             adr.GetHostName(), adr.GetHostAddress(), adr.GetPort());
1069       delete fSockets.RemoveAt(iSocket);
1070       continue;
1071     }
1072     if (strcmp(controlMessage, "ok") != 0) {
1073       TInetAddress adr = socket->GetInetAddress();
1074       Error("BroadcastHistos", "no \"ok\" message from client - "
1075             "disconnect client:\n %s (%s), port %d\n",
1076             adr.GetHostName(), adr.GetHostAddress(), adr.GetPort());
1077       delete fSockets.RemoveAt(iSocket);
1078       continue;
1079     }
1080
1081     socket->SetOption(kNoBlock, 0);
1082     if (socket->Send(message) < 0) {
1083       // remove the socket from the list if there was an error
1084       TInetAddress adr = socket->GetInetAddress();
1085       Error("BroadcastHistos", "sending histograms failed - "
1086             "disconnect client:\n %s (%s), port %d\n",
1087             adr.GetHostName(), adr.GetHostAddress(), adr.GetPort());
1088       delete fSockets.RemoveAt(iSocket);
1089     } else {
1090       gSystem->Sleep(100);
1091       socket->SetOption(kNoBlock, 1);
1092     }
1093   }
1094   fSockets.Compress();
1095 }
1096
1097
1098 //_____________________________________________________________________________
1099 AliMonitorProcess::AliMonitorInterruptHandler::AliMonitorInterruptHandler
1100   (AliMonitorProcess* process):
1101   TSignalHandler(kSigUser1, kFALSE), 
1102   fProcess(process) 
1103 {
1104 // constructor: set process
1105 }
1106
1107 //_____________________________________________________________________________
1108 AliMonitorProcess::AliMonitorInterruptHandler::AliMonitorInterruptHandler
1109   (const AliMonitorInterruptHandler& handler):
1110   TSignalHandler(handler)
1111 {
1112 // copy constructor
1113
1114   Fatal("AliMonitorInterruptHandler", "copy constructor not implemented");
1115 }
1116
1117 //_____________________________________________________________________________
1118 AliMonitorProcess::AliMonitorInterruptHandler& 
1119   AliMonitorProcess::AliMonitorInterruptHandler::operator = 
1120   (const AliMonitorInterruptHandler& /*handler*/) 
1121 {
1122 // assignment operator
1123
1124   Fatal("operator =", "assignment operator not implemented"); 
1125   return *this;
1126 }
1127
1128 //_____________________________________________________________________________
1129 Bool_t AliMonitorProcess::AliMonitorInterruptHandler::Notify() 
1130 {
1131 // interrupt signal -> stop process
1132
1133   Info("Notify", "the monitoring process will be stopped.");
1134   fProcess->Stop(); 
1135   return kTRUE;
1136 }