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