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