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