1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////////
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. //
26 ///////////////////////////////////////////////////////////////////////////////
30 #include <TGridResult.h>
33 #include <TServerSocket.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"
51 #include "AliTPCclustererMI.h"
52 #include "AliTPCtrackerMI.h"
53 #include "AliV0vertexer.h"
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>
70 ClassImp(AliMonitorProcess)
73 const Int_t AliMonitorProcess::fgkPort = 9327;
76 //_____________________________________________________________________________
77 AliMonitorProcess::AliMonitorProcess(
78 #if ROOT_VERSION_CODE <= 199169 // 3.10/01
79 const char* /*alienHost*/,
81 const char* alienHost,
84 const char* fileNameGalice)
86 // initialize the monitoring process and the monitor histograms
88 #if ROOT_VERSION_CODE <= 199169 // 3.10/01
89 fGrid = TGrid::Connect("alien", gSystem->Getenv("USER"));
91 fGrid = TGrid::Connect(alienHost, gSystem->Getenv("USER"));
93 if (!fGrid || fGrid->IsZombie() || !fGrid->IsConnected()) {
95 Fatal("AliMonitorProcess", "could not connect to alien");
97 #if ROOT_VERSION_CODE <= 199169 // 3.10/01
100 fAlienDir = alienDir;
102 fLogicalFileName = "";
105 fRunLoader = AliRunLoader::Open(fileNameGalice);
106 if (!fRunLoader) Fatal("AliMonitorProcess",
107 "could not get run loader from file %s",
110 fRunLoader->CdGAFile();
111 fTPCParam = AliTPC::LoadTPCParam(gFile);
112 if (!fTPCParam) Fatal("AliMonitorProcess", "could not load TPC parameters");
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");
123 // Init TPC parameters for HLT
124 Bool_t isinit=AliL3Transform::Init(const_cast<char*>(fileNameGalice),kTRUE);
126 cerr << "Could not create transform settings, please check log for error messages!" << endl;
135 fWriteHistoList = kFALSE;
137 fTopFolder = new TFolder("Monitor", "monitor histograms");
138 fTopFolder->SetOwner(kTRUE);
140 fMonitors.Add(new AliMonitorTPC(fTPCParam));
141 fMonitors.Add(new AliMonitorITS(fITSgeom));
142 fMonitors.Add(new AliMonitorV0s);
144 fMonitors.Add(new AliMonitorHLT(fTPCParam));
145 fMonitors.Add(new AliMonitorHLTHough(fTPCParam));
148 for (Int_t iMonitor = 0; iMonitor < fMonitors.GetEntriesFast(); iMonitor++) {
149 ((AliMonitor*) fMonitors[iMonitor])->CreateHistos(fTopFolder);
152 TIterator* iFolder = fTopFolder->GetListOfFolders()->MakeIterator();
153 while (TFolder* folder = (TFolder*) iFolder->Next()) folder->SetOwner(kTRUE);
156 fFile = TFile::Open("monitor_tree.root", "RECREATE");
157 if (!fFile || !fFile->IsOpen()) {
158 Fatal("AliMonitorProcess", "could not open file for tree");
160 fTree = new TTree("MonitorTree", "tree for monitoring");
161 for (Int_t iMonitor = 0; iMonitor < fMonitors.GetEntriesFast(); iMonitor++) {
162 ((AliMonitor*) fMonitors[iMonitor])->CreateBranches(fTree);
166 fServerSocket = new TServerSocket(fgkPort, kTRUE);
167 fServerSocket->SetOption(kNoBlock, 1);
168 fDisplaySocket = NULL;
169 CheckForConnections();
177 fInterruptHandler = new AliMonitorInterruptHandler(this);
178 gSystem->AddSignalHandler(fInterruptHandler);
181 //_____________________________________________________________________________
182 AliMonitorProcess::AliMonitorProcess(const AliMonitorProcess& process) :
185 Fatal("AliMonitorProcess", "copy constructor not implemented");
188 //_____________________________________________________________________________
189 AliMonitorProcess& AliMonitorProcess::operator = (const AliMonitorProcess&
192 Fatal("operator =", "assignment operator not implemented");
196 //_____________________________________________________________________________
197 AliMonitorProcess::~AliMonitorProcess()
205 delete fServerSocket;
207 delete fDisplaySocket;
209 #if ROOT_VERSION_CODE <= 199169 // 3.10/01
216 gSystem->Unlink("monitor_tree.root");
223 gSystem->RemoveSignalHandler(fInterruptHandler);
224 delete fInterruptHandler;
228 //_____________________________________________________________________________
229 const char* AliMonitorProcess::GetRevision()
235 //_____________________________________________________________________________
236 void AliMonitorProcess::SetStatus(EStatus status)
238 // set the current status and process system events
241 gSystem->ProcessEvents();
245 //_____________________________________________________________________________
246 void AliMonitorProcess::Run()
248 // run the monitor process:
249 // check for a raw data file, process the raw data file and delete it
255 while (!CheckForNewFile()) {
256 CheckForConnections();
258 if (fStopping) break;
261 if (fStopping) break;
273 //_____________________________________________________________________________
274 void AliMonitorProcess::Stop()
276 // set the fStopping flag to terminate the monitor process after the current
277 // event was processed
279 if (GetStatus() != kStopped) fStopping = kTRUE;
283 //_____________________________________________________________________________
284 void AliMonitorProcess::ProcessFile(const char* fileName)
286 // create a file with monitor histograms for a single file
288 if (GetStatus() != kStopped) {
289 Error("ProcessFile", "ProcessFile can not be called"
290 " while the monitor process is running");
294 fFileName = fileName;
295 Int_t nEventMin = fNEventsMin;
299 fNEventsMin = nEventMin;
304 //_____________________________________________________________________________
305 Bool_t AliMonitorProcess::CheckForNewFile()
307 // check whether a new file was registered in alien
309 #if ROOT_VERSION_CODE <= 199169 // 3.10/01
310 TGridResult* result = fGrid->Ls();
314 sprintf(dirName, "%s/adc-%d", fAlienDir.Data(), datime.GetDate());
316 sprintf(findName, "*.root");
317 Grid_ResultHandle_t handle = fGrid->Find(dirName, findName);
319 Error("CheckForNewFile", "could not open alien directory %s",
323 TGridResult* result = fGrid->CreateGridResult(handle);
329 #if ROOT_VERSION_CODE <= 199169 // 3.10/01
330 while (const char* entry = result->Next()) {
332 while (Grid_Result_t* resultEntry = result->Next()) {
333 const char* entry = resultEntry->name.c_str();
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);
347 if ((date > maxDate) || ((date == maxDate) && (time > maxTime))) {
355 if (maxDate < 0) return kFALSE; // no files found
356 if (fLogicalFileName.CompareTo(fileName) == 0) return kFALSE; // no new file
358 fLogicalFileName = fileName;
359 #if ROOT_VERSION_CODE <= 199169 // 3.10/01
360 result = fGrid->GetPhysicalFileNames(fLogicalFileName.Data());
361 fFileName = result->Next();
363 fileName = dirName + ("/" + fLogicalFileName);
364 handle = fGrid->GetPhysicalFileNames(fileName.Data());
366 Error("CheckForNewFile", "could not get physical file names for %s",
370 result = fGrid->CreateGridResult(handle);
372 Grid_Result_t* resultEntry = result->Next();
374 Error("CheckForNewFile", "could not get physical file names for %s",
378 fFileName = resultEntry->name2.c_str();
379 fFileName.ReplaceAll("castor:/", "rfio:/");
386 //_____________________________________________________________________________
387 Bool_t AliMonitorProcess::ProcessFile()
389 // loop over all events in the raw data file, run the reconstruction
390 // and fill the monitor histograms
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());
397 CreateHLT(fFileName);
398 CreateHLTHough(fFileName);
401 // loop over the events
402 for (Int_t iEvent = 0; iEvent < nEvents; iEvent++) {
403 CheckForConnections();
405 fRunLoader->SetEventNumber(0);
406 AliRawReaderRoot rawReader(fFileName, iEvent);
407 if (fStopping) break;
408 if (rawReader.GetRunNumber() != fRunNumber) {
411 fRunNumber = rawReader.GetRunNumber();
412 fEventNumber[0] = rawReader.GetEventId()[0];
413 fEventNumber[1] = rawReader.GetEventId()[1];
415 if (fStopping) break;
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]);
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;
441 if (fDisplaySocket) fDisplaySocket->Send("new event");
443 Info("ProcessFile", "filling histograms...");
444 for (Int_t iMonitor = 0; iMonitor < fMonitors.GetEntriesFast(); iMonitor++) {
445 CheckForConnections();
447 ((AliMonitor*) fMonitors[iMonitor])->FillHistos(fRunLoader, &rawReader,
449 if (fStopping) break;
451 if (fStopping) break;
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()) {
465 if (fStopping) break;
467 Info("ProcessFile", "filling the tree...");
470 Info("ProcessFile", "broadcasting histograms...");
471 CheckForConnections();
475 if (fStopping) break;
486 //_____________________________________________________________________________
487 void AliMonitorProcess::Reset()
489 // write the current histograms to a file and reset them
491 if (fSubRunNumber == 0) fSubRunNumber++;
498 //_____________________________________________________________________________
499 UInt_t AliMonitorProcess::GetEventPeriodNumber() const
501 // get the period number from the event id
503 return (fEventNumber[1] >> 4);
506 //_____________________________________________________________________________
507 UInt_t AliMonitorProcess::GetEventOrbitNumber() const
509 // get the orbit number from the event id
511 return ((fEventNumber[1] & 0x000F) << 20) + (fEventNumber[0] >> 12);
514 //_____________________________________________________________________________
515 UInt_t AliMonitorProcess::GetEventBunchNumber() const
517 // get the bunch number from the event id
519 return (fEventNumber[0] % 0x0FFF);
522 //_____________________________________________________________________________
523 Int_t AliMonitorProcess::GetNumberOfEvents(const char* fileName) const
525 // determine the number of events in the given raw data file
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;
536 TTree* tree = (TTree*) file->Get("RAW");
538 Error("GetNumberOfEvents", "could not find tree with raw data");
540 nEvents = (Int_t) tree->GetEntries();
548 //_____________________________________________________________________________
549 Bool_t AliMonitorProcess::ReconstructTPC(AliRawReader* rawReader, AliESD* esd)
551 // find TPC clusters and tracks
555 AliLoader* tpcLoader = fRunLoader->GetLoader("TPCLoader");
557 Error("ReconstructTPC", "no TPC loader found");
560 gSystem->Unlink("TPC.RecPoints.root");
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");
572 Info("ReconstructTPC", "reconstructing tracks...");
573 AliTPCtrackerMI tracker(fTPCParam);
574 tracker.LoadClusters(tpcLoader->TreeR());
575 tracker.Clusters2Tracks(esd);
576 tracker.UnloadClusters();
577 tpcLoader->UnloadRecPoints();
582 //_____________________________________________________________________________
583 Bool_t AliMonitorProcess::ReconstructITS(AliRawReader* rawReader, AliESD* esd)
585 // find ITS clusters and tracks
589 AliLoader* itsLoader = fRunLoader->GetLoader("ITSLoader");
591 Error("ReconstructITS", "no ITS loader found");
594 gSystem->Unlink("ITS.RecPoints.root");
597 Info("ReconstructITS", "reconstructing clusters...");
598 itsLoader->LoadRecPoints("recreate");
599 AliITSclustererV2 clusterer(fITSgeom);
600 itsLoader->MakeRecPointsContainer();
601 clusterer.Digits2Clusters(rawReader);
604 Info("ReconstructITS", "reconstructing tracks...");
605 AliITStrackerV2 tracker(fITSgeom);
606 tracker.LoadClusters(itsLoader->TreeR());
607 tracker.Clusters2Tracks(esd);
608 tracker.UnloadClusters();
610 itsLoader->UnloadRecPoints();
614 //_____________________________________________________________________________
615 Bool_t AliMonitorProcess::ReconstructV0s(AliESD* esd)
622 Info("ReconstructV0s", "reconstructing V0s...");
623 AliV0vertexer vertexer;
625 esd->GetVertex()->GetXYZ(vtx);
626 vertexer.SetVertex(vtx);
627 vertexer.Tracks2V0vertices(esd);
632 //_____________________________________________________________________________
634 void AliMonitorProcess::CreateHLT(const char* fileName)
637 // create the HLT (Level3) object
639 if (fHLT) delete fHLT;
642 strcpy(name, fileName);
643 fHLT = new AliLevel3(name);
644 fHLT->Init("./", AliLevel3::kRaw, 1);
646 fHLT->SetClusterFinderParam(-1, -1, kTRUE);
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);
668 fHLT->WriteFiles("./hlt/");
671 //_____________________________________________________________________________
672 void AliMonitorProcess::CreateHLTHough(const char* fileName)
675 // create the HLT Hough transform (L3Hough) object
677 if (fHLTHough) delete fHLTHough;
680 strcpy(name, fileName);
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);
693 //_____________________________________________________________________________
694 Bool_t AliMonitorProcess::ReconstructHLT(
702 // run the HLT cluster and track finder
707 Warning("ReconstructHLT", "the code was compiled without HLT support");
711 gSystem->Exec("rm -rf hlt");
712 gSystem->MakeDirectory("hlt");
713 if (!fHLT) return kFALSE;
715 fHLT->ProcessEvent(0, 35, iEvent);
717 // remove the event number from the file names
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);
731 //_____________________________________________________________________________
732 Bool_t AliMonitorProcess::ReconstructHLTHough(
740 // run the HLT Hough transformer
745 Warning("ReconstructHLTHough", "the code was compiled without HLT support");
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;
755 // fHLTHough->Process(0, 35);
756 // Loop over TPC sectors and process the data
757 for(Int_t i=0; i<=35; i++)
759 fHLTHough->ReadData(i,iEvent);
760 fHLTHough->Transform();
761 // if(fHLTHough->fAddHistograms)
762 fHLTHough->AddAllHistograms();
763 fHLTHough->FindTrackCandidates();
764 fHLTHough->AddTracks();
766 fHLTHough->WriteTracks("./hlt/hough");
768 // Run cluster fitter
769 AliL3ClusterFitter *fitter = new AliL3ClusterFitter("./hlt");
771 // Set debug flag for the cluster fitter
774 // Setting fitter parameters
775 fitter->SetInnerWidthFactor(1,1.5);
776 fitter->SetOuterWidthFactor(1,1.5);
777 fitter->SetNmaxOverlaps(5);
779 //fitter->SetChiSqMax(5,kFALSE); //isolated clusters
780 fitter->SetChiSqMax(5,kTRUE); //overlapping clusters
782 Int_t rowrange[2] = {0,AliL3Transform::GetNRows()-1};
784 // Takes input from global hough tracks produced by HT
785 fitter->LoadSeeds(rowrange,kFALSE,iEvent);
789 for(Int_t islice = 0; islice <= 35; islice++)
791 for(Int_t ipatch = 0; ipatch < AliL3Transform::GetNPatches(); ipatch++)
794 fHLTHough->GetMemHandler(ipatch)->Free();
795 fHLTHough->GetMemHandler(ipatch)->Init(islice,ipatch);
796 AliL3DigitRowData *digits = (AliL3DigitRowData *)fHLTHough->GetMemHandler(ipatch)->AliAltroDigits2Memory(ndigits,iEvent);
798 fitter->Init(islice,ipatch);
799 fitter->SetInputData(digits);
800 fitter->FindClusters();
801 fitter->WriteClusters();
805 // Refit of the clusters
807 //The seeds are the input tracks from circle HT
808 AliL3TrackArray *tracks = fitter->GetSeeds();
809 AliL3Fitter *ft = new AliL3Fitter(&vertex,1);
811 ft->LoadClusters("./hlt/fitter/",iEvent,kFALSE);
812 for(Int_t i=0; i<tracks->GetNTracks(); i++)
814 AliL3Track *track = tracks->GetCheckedTrack(i);
816 if(track->GetNHits() < 20) continue;
817 ft->SortTrackClusters(track);
819 track->UpdateToFirstPoint();
823 //Write the final tracks
824 fitter->WriteTracks(20);
828 // remove the event number from the file names
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);
840 //_____________________________________________________________________________
841 Bool_t AliMonitorProcess::WriteHistos()
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
849 // rename tree file and create a new one
856 sprintf(fileName, "monitor_tree_%d.root", fRunNumber);
857 if (fSubRunNumber > 0) {
858 sprintf(fileName, "monitor_tree_%d_%d.root", fRunNumber, fSubRunNumber);
860 if (fNEvents < fNEventsMin) {
861 gSystem->Unlink("monitor_tree.root");
863 gSystem->Rename("monitor_tree.root", fileName);
866 fFile = TFile::Open("monitor_tree.root", "RECREATE");
867 if (!fFile || !fFile->IsOpen()) {
868 Fatal("WriteHistos", "could not open file for tree");
870 fTree = new TTree("MonitorTree", "tree for monitoring");
871 for (Int_t iMonitor = 0; iMonitor < fMonitors.GetEntriesFast(); iMonitor++) {
872 ((AliMonitor*) fMonitors[iMonitor])->CreateBranches(fTree);
876 // write the histograms
877 if (fNEvents < fNEventsMin) return kTRUE;
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()) {
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);
896 TFile* file = TFile::Open(fileName, "recreate");
897 if (!file || !file->IsOpen()) {
898 Error("WriteHistos", "could not open file %s", fileName);
904 if (file) delete file;
909 //_____________________________________________________________________________
910 void AliMonitorProcess::StartNewRun()
912 // reset the histograms for a new run
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()) {
929 //_____________________________________________________________________________
930 void AliMonitorProcess::CheckForConnections()
932 // check if new clients want to connect and add them to the list of sockets
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());
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;
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());
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());
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);
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);
1001 //_____________________________________________________________________________
1002 void AliMonitorProcess::BroadcastHistos(TSocket* toSocket)
1004 // send the monitor histograms to the clients
1006 SetStatus(kBroadcasting);
1007 TMessage message(kMESS_OBJECT);
1008 message.WriteObject(fTopFolder);
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;
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);
1024 // receive control message
1025 char controlMessage[256];
1026 Int_t result = socket->Recv(controlMessage, 255);
1028 gSystem->Sleep(1000); // wait one second and try again
1029 result = socket->Recv(controlMessage, 255);
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);
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);
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);
1057 gSystem->Sleep(100);
1058 socket->SetOption(kNoBlock, 1);
1061 fSockets.Compress();
1065 //_____________________________________________________________________________
1066 AliMonitorProcess::AliMonitorInterruptHandler::AliMonitorInterruptHandler
1067 (AliMonitorProcess* process):
1068 TSignalHandler(kSigUser1, kFALSE),
1071 // constructor: set process
1074 //_____________________________________________________________________________
1075 AliMonitorProcess::AliMonitorInterruptHandler::AliMonitorInterruptHandler
1076 (const AliMonitorInterruptHandler& handler):
1077 TSignalHandler(handler)
1081 Fatal("AliMonitorInterruptHandler", "copy constructor not implemented");
1084 //_____________________________________________________________________________
1085 AliMonitorProcess::AliMonitorInterruptHandler&
1086 AliMonitorProcess::AliMonitorInterruptHandler::operator =
1087 (const AliMonitorInterruptHandler& /*handler*/)
1089 // assignment operator
1091 Fatal("operator =", "assignment operator not implemented");
1095 //_____________________________________________________________________________
1096 Bool_t AliMonitorProcess::AliMonitorInterruptHandler::Notify()
1098 // interrupt signal -> stop process
1100 Info("Notify", "the monitoring process will be stopped.");