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