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 ///////////////////////////////////////////////////////////////////////////////
29 #include "AliMonitorProcess.h"
30 #include "AliMonitorTPC.h"
31 #include "AliMonitorITS.h"
32 #include "AliMonitorV0s.h"
33 #include "AliMonitorHLT.h"
34 #include "AliMonitorHLTHough.h"
35 #include "AliRawReaderRoot.h"
36 #include "AliLoader.h"
39 #include "AliTPCclustererMI.h"
40 #include "AliTPCtrackerMI.h"
42 #include "AliITSclustererV2.h"
43 #include "AliITStrackerV2.h"
44 #include "AliITSLoader.h"
45 #include "AliV0vertexer.h"
48 #include <TServerSocket.h>
50 #include <TGridResult.h>
53 #include <AliLevel3.h>
54 #include <AliL3Transform.h>
55 #include <AliL3Track.h>
56 #include <AliL3TrackArray.h>
57 #include <AliL3StandardIncludes.h>
58 #include <AliL3HoughMaxFinder.h>
59 #include <AliL3HoughBaseTransformer.h>
60 #include <AliL3Hough.h>
61 #include <AliL3ClusterFitter.h>
62 #include <AliL3Vertex.h>
63 #include <AliL3Fitter.h>
64 #include <AliL3DDLDataFileHandler.h>
67 ClassImp(AliMonitorProcess)
70 const Int_t AliMonitorProcess::fgkPort = 9327;
73 //_____________________________________________________________________________
74 AliMonitorProcess::AliMonitorProcess(
75 #if ROOT_VERSION_CODE <= 199169 // 3.10/01
76 const char* /*alienHost*/,
78 const char* alienHost,
81 const char* fileNameGalice)
83 // initialize the monitoring process and the monitor histograms
85 #if ROOT_VERSION_CODE <= 199169 // 3.10/01
86 fGrid = TGrid::Connect("alien", gSystem->Getenv("USER"));
88 fGrid = TGrid::Connect(alienHost, gSystem->Getenv("USER"));
90 if (!fGrid || fGrid->IsZombie() || !fGrid->IsConnected()) {
92 Fatal("AliMonitorProcess", "could not connect to alien");
94 #if ROOT_VERSION_CODE <= 199169 // 3.10/01
99 fLogicalFileName = "";
102 fRunLoader = AliRunLoader::Open(fileNameGalice);
103 if (!fRunLoader) Fatal("AliMonitorProcess",
104 "could not get run loader from file %s",
107 fRunLoader->CdGAFile();
108 fTPCParam = AliTPC::LoadTPCParam(gFile);
109 if (!fTPCParam) Fatal("AliMonitorProcess", "could not load TPC parameters");
111 fRunLoader->LoadgAlice();
112 gAlice = fRunLoader->GetAliRun();
113 if (!gAlice) Fatal("AliMonitorProcess", "no gAlice object found");
114 AliITS* its = (AliITS*) gAlice->GetModule("ITS");
115 if (!its) Fatal("AliMonitorProcess", "no ITS detector found");
116 fITSgeom = its->GetITSgeom();
117 if (!fITSgeom) Fatal("AliMonitorProcess", "could not load ITS geometry");
120 // Init TPC parameters for HLT
121 Bool_t isinit=AliL3Transform::Init(const_cast<char*>(fileNameGalice),kTRUE);
123 cerr << "Could not create transform settings, please check log for error messages!" << endl;
132 fWriteHistoList = kFALSE;
134 fTopFolder = new TFolder("Monitor", "monitor histograms");
135 fTopFolder->SetOwner(kTRUE);
137 fMonitors.Add(new AliMonitorTPC(fTPCParam));
138 fMonitors.Add(new AliMonitorITS(fITSgeom));
139 fMonitors.Add(new AliMonitorV0s);
141 fMonitors.Add(new AliMonitorHLT(fTPCParam));
142 fMonitors.Add(new AliMonitorHLTHough(fTPCParam));
145 for (Int_t iMonitor = 0; iMonitor < fMonitors.GetEntriesFast(); iMonitor++) {
146 ((AliMonitor*) fMonitors[iMonitor])->CreateHistos(fTopFolder);
149 TIterator* iFolder = fTopFolder->GetListOfFolders()->MakeIterator();
150 while (TFolder* folder = (TFolder*) iFolder->Next()) folder->SetOwner(kTRUE);
153 fFile = TFile::Open("monitor_tree.root", "RECREATE");
154 if (!fFile || !fFile->IsOpen()) {
155 Fatal("AliMonitorProcess", "could not open file for tree");
157 fTree = new TTree("MonitorTree", "tree for monitoring");
158 for (Int_t iMonitor = 0; iMonitor < fMonitors.GetEntriesFast(); iMonitor++) {
159 ((AliMonitor*) fMonitors[iMonitor])->CreateBranches(fTree);
163 fServerSocket = new TServerSocket(fgkPort, kTRUE);
164 fServerSocket->SetOption(kNoBlock, 1);
165 fDisplaySocket = NULL;
166 CheckForConnections();
175 //_____________________________________________________________________________
176 AliMonitorProcess::AliMonitorProcess(const AliMonitorProcess& process) :
179 Fatal("AliMonitorProcess", "copy constructor not implemented");
182 //_____________________________________________________________________________
183 AliMonitorProcess& AliMonitorProcess::operator = (const AliMonitorProcess&
186 Fatal("operator =", "assignment operator not implemented");
190 //_____________________________________________________________________________
191 AliMonitorProcess::~AliMonitorProcess()
199 delete fServerSocket;
201 delete fDisplaySocket;
203 #if ROOT_VERSION_CODE <= 199169 // 3.10/01
210 gSystem->Unlink("monitor_tree.root");
219 //_____________________________________________________________________________
220 const char* AliMonitorProcess::GetRevision()
226 //_____________________________________________________________________________
227 void AliMonitorProcess::SetStatus(EStatus status)
229 // set the current status and process system events
232 gSystem->ProcessEvents();
236 //_____________________________________________________________________________
237 void AliMonitorProcess::Run()
239 // run the monitor process:
240 // check for a raw data file, process the raw data file and delete it
246 while (!CheckForNewFile()) {
247 CheckForConnections();
249 if (fStopping) break;
252 if (fStopping) break;
264 //_____________________________________________________________________________
265 void AliMonitorProcess::Stop()
267 // set the fStopping flag to terminate the monitor process after the current
268 // event was processed
270 if (GetStatus() != kStopped) fStopping = kTRUE;
274 //_____________________________________________________________________________
275 void AliMonitorProcess::ProcessFile(const char* fileName)
277 // create a file with monitor histograms for a single file
279 if (GetStatus() != kStopped) {
280 Error("ProcessFile", "ProcessFile can not be called"
281 " while the monitor process is running");
285 fFileName = fileName;
286 Int_t nEventMin = fNEventsMin;
290 fNEventsMin = nEventMin;
295 //_____________________________________________________________________________
296 Bool_t AliMonitorProcess::CheckForNewFile()
298 // check whether a new file was registered in alien
300 #if ROOT_VERSION_CODE <= 199169 // 3.10/01
301 TGridResult* result = fGrid->Ls();
303 Grid_ResultHandle_t handle = fGrid->OpenDir(fAlienDir);
305 Error("CheckForNewFile", "could not open alien directory %s",
309 TGridResult* result = fGrid->CreateGridResult(handle);
315 #if ROOT_VERSION_CODE <= 199169 // 3.10/01
316 while (const char* entry = result->Next()) {
318 while (Grid_Result_t* resultEntry = result->Next()) {
319 const char* entry = resultEntry->name.c_str();
321 // entry = host_date_time.root
322 TString entryCopy(entry);
323 char* p = const_cast<char*>(entryCopy.Data());
324 if (!strtok(p, "_") || !p) continue; // host name
325 char* dateStr = strtok(NULL, "_");
326 if (!dateStr || !p) continue;
327 char* timeStr = strtok(NULL, ".");
328 if (!timeStr || !p) continue;
329 Long_t date = atoi(dateStr);
330 Long_t time = atoi(timeStr);
332 if ((date > maxDate) || ((date == maxDate) && (time > maxTime))) {
340 if (maxDate < 0) return kFALSE; // no files found
341 if (fLogicalFileName.CompareTo(fileName) == 0) return kFALSE; // no new file
343 fLogicalFileName = fileName;
344 #if ROOT_VERSION_CODE <= 199169 // 3.10/01
345 result = fGrid->GetPhysicalFileNames(fLogicalFileName.Data());
346 fFileName = result->Next();
348 fileName = fAlienDir + "/" + fLogicalFileName;
349 handle = fGrid->GetPhysicalFileNames(fileName.Data());
351 Error("CheckForNewFile", "could not get physical file names for %s",
355 result = fGrid->CreateGridResult(handle);
357 Grid_Result_t* resultEntry = result->Next();
359 Error("CheckForNewFile", "could not get physical file names for %s",
363 fFileName = resultEntry->name2.c_str();
370 //_____________________________________________________________________________
371 Bool_t AliMonitorProcess::ProcessFile()
373 // loop over all events in the raw data file, run the reconstruction
374 // and fill the monitor histograms
376 Int_t nEvents = GetNumberOfEvents(fFileName);
377 if (nEvents <= 0) return kFALSE;
378 Info("ProcessFile", "found %d event(s) in file %s",
379 nEvents, fFileName.Data());
381 CreateHLT(fFileName);
382 CreateHLTHough(fFileName);
385 // loop over the events
386 for (Int_t iEvent = 0; iEvent < nEvents; iEvent++) {
388 fRunLoader->SetEventNumber(0);
389 AliRawReaderRoot rawReader(fFileName, iEvent);
390 if (fStopping) break;
391 if (rawReader.GetRunNumber() != fRunNumber) {
394 fRunNumber = rawReader.GetRunNumber();
395 fEventNumber[0] = rawReader.GetEventId()[0];
396 fEventNumber[1] = rawReader.GetEventId()[1];
398 if (fStopping) break;
401 if (!ReconstructTPC(&rawReader)) return kFALSE;
402 if (fStopping) break;
403 if (!ReconstructITS(&rawReader)) return kFALSE;
404 if (fStopping) break;
405 if (!ReconstructV0s()) return kFALSE;
406 if (fStopping) break;
407 if (!ReconstructHLT(iEvent)) return kFALSE;
408 if (fStopping) break;
409 if (!ReconstructHLTHough(iEvent)) return kFALSE;
410 if (fStopping) break;
412 if (fDisplaySocket) fDisplaySocket->Send("new event");
414 Info("ProcessFile", "filling histograms...");
416 for (Int_t iMonitor = 0; iMonitor < fMonitors.GetEntriesFast(); iMonitor++) {
417 ((AliMonitor*) fMonitors[iMonitor])->FillHistos(fRunLoader, &rawReader);
418 if (fStopping) break;
420 if (fStopping) break;
422 Info("ProcessFile", "updating histograms...");
423 SetStatus(kUpdating);
424 TIterator* iFolder = fTopFolder->GetListOfFolders()->MakeIterator();
425 while (TFolder* folder = (TFolder*) iFolder->Next()) {
426 TIterator* iHisto = folder->GetListOfFolders()->MakeIterator();
427 while (AliMonitorPlot* histo = (AliMonitorPlot*) iHisto->Next()) {
433 if (fStopping) break;
435 Info("ProcessFile", "filling the tree...");
438 Info("ProcessFile", "broadcasting histograms...");
439 CheckForConnections();
443 if (fStopping) break;
454 //_____________________________________________________________________________
455 void AliMonitorProcess::Reset()
457 // write the current histograms to a file and reset them
459 if (fSubRunNumber == 0) fSubRunNumber++;
466 //_____________________________________________________________________________
467 UInt_t AliMonitorProcess::GetEventPeriodNumber() const
469 // get the period number from the event id
471 return (fEventNumber[1] >> 4);
474 //_____________________________________________________________________________
475 UInt_t AliMonitorProcess::GetEventOrbitNumber() const
477 // get the orbit number from the event id
479 return ((fEventNumber[1] & 0x000F) << 20) + (fEventNumber[0] >> 12);
482 //_____________________________________________________________________________
483 UInt_t AliMonitorProcess::GetEventBunchNumber() const
485 // get the bunch number from the event id
487 return (fEventNumber[0] % 0x0FFF);
490 //_____________________________________________________________________________
491 Int_t AliMonitorProcess::GetNumberOfEvents(const char* fileName) const
493 // determine the number of events in the given raw data file
497 TFile* file = TFile::Open(fileName);
498 if (!file || !file->IsOpen()) {
499 Error("GetNumberOfEvents", "could not open file %s", fileName);
500 if (file) delete file;
504 TTree* tree = (TTree*) file->Get("RAW");
506 Error("GetNumberOfEvents", "could not find tree with raw data");
508 nEvents = (Int_t) tree->GetEntries();
516 //_____________________________________________________________________________
517 Bool_t AliMonitorProcess::ReconstructTPC(AliRawReader* rawReader)
519 // find TPC clusters and tracks
523 AliLoader* tpcLoader = fRunLoader->GetLoader("TPCLoader");
525 Error("ReconstructTPC", "no TPC loader found");
528 gSystem->Unlink("TPC.RecPoints.root");
529 gSystem->Unlink("TPC.Tracks.root");
532 Info("ReconstructTPC", "reconstructing clusters...");
533 tpcLoader->LoadRecPoints("recreate");
534 AliTPCclustererMI clusterer(fTPCParam);
535 tpcLoader->MakeRecPointsContainer();
536 clusterer.SetOutput(tpcLoader->TreeR());
537 clusterer.Digits2Clusters(rawReader);
538 tpcLoader->WriteRecPoints("OVERWRITE");
541 Info("ReconstructTPC", "reconstructing tracks...");
542 tpcLoader->LoadTracks("recreate");
544 AliTPCtrackerMI tracker(fTPCParam);
545 tracker.Clusters2Tracks();
548 tpcLoader->UnloadRecPoints();
549 tpcLoader->UnloadTracks();
553 //_____________________________________________________________________________
554 Bool_t AliMonitorProcess::ReconstructITS(AliRawReader* rawReader)
556 // find ITS clusters and tracks
560 AliLoader* itsLoader = fRunLoader->GetLoader("ITSLoader");
562 Error("ReconstructITS", "no ITS loader found");
565 AliLoader* tpcLoader = fRunLoader->GetLoader("TPCLoader");
567 Error("ReconstructITS", "no TPC loader found");
570 gSystem->Unlink("ITS.RecPoints.root");
571 gSystem->Unlink("ITS.Tracks.root");
574 Info("ReconstructITS", "reconstructing clusters...");
575 itsLoader->LoadRecPoints("recreate");
576 AliITSclustererV2 clusterer(fITSgeom);
577 itsLoader->MakeRecPointsContainer();
578 clusterer.Digits2Clusters(rawReader);
581 Info("ReconstructITS", "reconstructing tracks...");
582 itsLoader->LoadTracks("recreate");
583 itsLoader->MakeTracksContainer();
584 tpcLoader->LoadTracks();
585 AliITStrackerV2 tracker(fITSgeom);
586 tracker.LoadClusters(itsLoader->TreeR());
587 tracker.Clusters2Tracks(tpcLoader->TreeT(), itsLoader->TreeT());
588 tracker.UnloadClusters();
589 itsLoader->WriteTracks("OVERWRITE");
591 itsLoader->UnloadRecPoints();
592 itsLoader->UnloadTracks();
593 tpcLoader->UnloadTracks();
597 //_____________________________________________________________________________
598 Bool_t AliMonitorProcess::ReconstructV0s()
604 AliITSLoader* itsLoader = (AliITSLoader*) fRunLoader->GetLoader("ITSLoader");
606 Error("ReconstructV0", "no ITS loader found");
609 gSystem->Unlink("ITS.V0s.root");
612 Info("ReconstructV0s", "reconstructing V0s...");
613 itsLoader->LoadTracks("read");
614 itsLoader->LoadV0s("recreate");
615 AliV0vertexer vertexer;
616 TTree* tracks = itsLoader->TreeT();
618 Error("ReconstructV0s", "no ITS tracks tree found");
621 if (!itsLoader->TreeV0()) itsLoader->MakeTree("V0");
622 TTree* v0s = itsLoader->TreeV0();
623 vertexer.Tracks2V0vertices(tracks, v0s);
624 itsLoader->WriteV0s("OVERWRITE");
626 itsLoader->UnloadTracks();
627 itsLoader->UnloadV0s();
631 //_____________________________________________________________________________
633 void AliMonitorProcess::CreateHLT(const char* fileName)
636 // create the HLT (Level3) object
638 if (fHLT) delete fHLT;
641 strcpy(name, fileName);
642 fHLT = new AliLevel3(name);
643 fHLT->Init("./", AliLevel3::kRaw, 1);
645 fHLT->SetClusterFinderParam(-1, -1, kTRUE);
647 Int_t phiSegments = 50;
648 Int_t etaSegments = 100;
649 Int_t trackletlength = 3;
650 Int_t tracklength = 20;//40 or 5
651 Int_t rowscopetracklet = 2;
652 Int_t rowscopetrack = 10;
653 Double_t minPtFit = 0;
654 Double_t maxangle = 0.1745;
655 Double_t goodDist = 5;
656 Double_t maxphi = 0.1;
657 Double_t maxeta = 0.1;
658 Double_t hitChi2Cut = 15;//100 or 15
659 Double_t goodHitChi2 = 5;//20 or 5
660 Double_t trackChi2Cut = 10;//50 or 10
661 fHLT->SetTrackerParam(phiSegments, etaSegments,
662 trackletlength, tracklength,
663 rowscopetracklet, rowscopetrack,
664 minPtFit, maxangle, goodDist, hitChi2Cut,
665 goodHitChi2, trackChi2Cut, 50, maxphi, maxeta, kTRUE);
667 fHLT->WriteFiles("./hlt/");
670 //_____________________________________________________________________________
671 void AliMonitorProcess::CreateHLTHough(const char* fileName)
674 // create the HLT Hough transform (L3Hough) object
676 if (fHLTHough) delete fHLTHough;
679 strcpy(name, fileName);
681 fHLTHough = new AliL3Hough();
682 fHLTHough->SetThreshold(4);
683 fHLTHough->SetTransformerParams(140,150,0.5,-1);
684 fHLTHough->SetPeakThreshold(9000,-1);// or 6000
685 fHLTHough->Init("./", kFALSE, 50, kFALSE,0,name);
686 fHLTHough->SetAddHistograms();
687 // fHLTHough->GetMaxFinder()->SetThreshold(14000);
692 //_____________________________________________________________________________
693 Bool_t AliMonitorProcess::ReconstructHLT(
701 // run the HLT cluster and track finder
706 Warning("ReconstructHLT", "the code was compiled without HLT support");
710 gSystem->Exec("rm -rf hlt");
711 gSystem->MakeDirectory("hlt");
712 if (!fHLT) return kFALSE;
714 fHLT->ProcessEvent(0, 35, iEvent);
716 // remove the event number from the file names
718 sprintf(command, "rename points_%d points hlt/*.raw", iEvent);
719 gSystem->Exec(command);
720 sprintf(command, "rename tracks_tr_%d tracks_tr hlt/*.raw", iEvent);
721 gSystem->Exec(command);
722 sprintf(command, "rename tracks_gl_%d tracks_gl hlt/*.raw", iEvent);
723 gSystem->Exec(command);
724 sprintf(command, "rename tracks_%d tracks hlt/*.raw", iEvent);
725 gSystem->Exec(command);
730 //_____________________________________________________________________________
731 Bool_t AliMonitorProcess::ReconstructHLTHough(
739 // run the HLT Hough transformer
744 Warning("ReconstructHLTHough", "the code was compiled without HLT support");
748 gSystem->Exec("rm -rf hlt/hough");
749 gSystem->MakeDirectory("hlt/hough");
750 gSystem->Exec("rm -rf hlt/fitter");
751 gSystem->MakeDirectory("hlt/fitter");
752 if (!fHLTHough) return kFALSE;
754 // fHLTHough->Process(0, 35);
755 // Loop over TPC sectors and process the data
756 for(Int_t i=0; i<=35; i++)
758 fHLTHough->ReadData(i,iEvent);
759 fHLTHough->Transform();
760 // if(fHLTHough->fAddHistograms)
761 fHLTHough->AddAllHistograms();
762 fHLTHough->FindTrackCandidates();
763 fHLTHough->AddTracks();
765 fHLTHough->WriteTracks("./hlt/hough");
767 // Run cluster fitter
768 AliL3ClusterFitter *fitter = new AliL3ClusterFitter("./hlt");
770 // Set debug flag for the cluster fitter
773 // Setting fitter parameters
774 fitter->SetInnerWidthFactor(1,1.5);
775 fitter->SetOuterWidthFactor(1,1.5);
776 fitter->SetNmaxOverlaps(5);
778 //fitter->SetChiSqMax(5,kFALSE); //isolated clusters
779 fitter->SetChiSqMax(5,kTRUE); //overlapping clusters
781 Int_t rowrange[2] = {0,AliL3Transform::GetNRows()-1};
783 // Takes input from global hough tracks produced by HT
784 fitter->LoadSeeds(rowrange,kFALSE,iEvent);
788 for(Int_t islice = 0; islice <= 35; islice++)
790 for(Int_t ipatch = 0; ipatch < AliL3Transform::GetNPatches(); ipatch++)
793 fHLTHough->GetMemHandler(ipatch)->Init(islice,ipatch);
794 AliL3DigitRowData *digits = (AliL3DigitRowData *)fHLTHough->GetMemHandler(ipatch)->AliAltroDigits2Memory(ndigits,iEvent);
796 fitter->Init(islice,ipatch);
797 fitter->SetInputData(digits);
798 fitter->FindClusters();
799 fitter->WriteClusters();
803 // Refit of the clusters
805 //The seeds are the input tracks from circle HT
806 AliL3TrackArray *tracks = fitter->GetSeeds();
807 AliL3Fitter *ft = new AliL3Fitter(&vertex,1);
809 ft->LoadClusters("./hlt/fitter/",iEvent,kFALSE);
810 for(Int_t i=0; i<tracks->GetNTracks(); i++)
812 AliL3Track *track = tracks->GetCheckedTrack(i);
814 if(track->GetNHits() < 20) continue;
815 ft->SortTrackClusters(track);
817 ft->UpdateTrack(track);
821 //Write the final tracks
822 fitter->WriteTracks(20);
826 // remove the event number from the file names
828 sprintf(command, "rename tracks_%d tracks hlt/hough/*.raw", iEvent);
829 gSystem->Exec(command);
830 sprintf(command, "rename tracks_%d tracks hlt/fitter/*.raw", iEvent);
831 gSystem->Exec(command);
832 sprintf(command, "rename points_%d points hlt/fitter/*.raw", iEvent);
833 gSystem->Exec(command);
838 //_____________________________________________________________________________
839 Bool_t AliMonitorProcess::WriteHistos()
841 // write the monitor tree and the monitor histograms to the file
842 // "monitor_<run number>[_<sub_run_number>].root"
843 // if at least fNEventsMin events were monitored
847 // rename tree file and create a new one
854 sprintf(fileName, "monitor_tree_%d.root", fRunNumber);
855 if (fSubRunNumber > 0) {
856 sprintf(fileName, "monitor_tree_%d_%d.root", fRunNumber, fSubRunNumber);
858 if (fNEvents < fNEventsMin) {
859 gSystem->Unlink("monitor_tree.root");
861 gSystem->Rename("monitor_tree.root", fileName);
864 fFile = TFile::Open("monitor_tree.root", "RECREATE");
865 if (!fFile || !fFile->IsOpen()) {
866 Fatal("WriteHistos", "could not open file for tree");
868 fTree = new TTree("MonitorTree", "tree for monitoring");
869 for (Int_t iMonitor = 0; iMonitor < fMonitors.GetEntriesFast(); iMonitor++) {
870 ((AliMonitor*) fMonitors[iMonitor])->CreateBranches(fTree);
874 // write the histograms
875 if (fNEvents < fNEventsMin) return kTRUE;
877 if (!fWriteHistoList) {
878 TIterator* iFolder = fTopFolder->GetListOfFolders()->MakeIterator();
879 while (TFolder* folder = (TFolder*) iFolder->Next()) {
880 TIterator* iHisto = folder->GetListOfFolders()->MakeIterator();
881 while (AliMonitorPlot* histo = (AliMonitorPlot*) iHisto->Next()) {
889 Bool_t result = kTRUE;
890 sprintf(fileName, "monitor_%d.root", fRunNumber);
891 if (fSubRunNumber > 0) {
892 sprintf(fileName, "monitor_%d_%d.root", fRunNumber, fSubRunNumber);
894 TFile* file = TFile::Open(fileName, "recreate");
895 if (!file || !file->IsOpen()) {
896 Error("WriteHistos", "could not open file %s", fileName);
902 if (file) delete file;
907 //_____________________________________________________________________________
908 void AliMonitorProcess::StartNewRun()
910 // reset the histograms for a new run
912 SetStatus(kResetting);
913 TIterator* iFolder = fTopFolder->GetListOfFolders()->MakeIterator();
914 while (TFolder* folder = (TFolder*) iFolder->Next()) {
915 TIterator* iHisto = folder->GetListOfFolders()->MakeIterator();
916 while (AliMonitorPlot* histo = (AliMonitorPlot*) iHisto->Next()) {
927 //_____________________________________________________________________________
928 void AliMonitorProcess::CheckForConnections()
930 // check if new clients want to connect and add them to the list of sockets
932 TMessage message(kMESS_OBJECT);
933 message.WriteObject(fTopFolder);
934 SetStatus(kConnecting);
937 while ((socket = fServerSocket->Accept()) != (TSocket*)-1) {
938 socket->SetOption(kNoBlock, 1);
939 char socketType[256];
940 if (!socket->Recv(socketType, 255)) continue;
941 if (strcmp(socketType, "client") == 0) {
942 if ((fNEvents == 0) || (socket->Send(message) >= 0)) {
943 fSockets.Add(socket);
944 TInetAddress adr = socket->GetInetAddress();
945 Info("CheckForConnections", "new client:\n %s (%s), port %d\n",
946 adr.GetHostName(), adr.GetHostAddress(), adr.GetPort());
948 } else if (strcmp(socketType, "display") == 0) {
949 if (fDisplaySocket) {
950 fDisplaySocket->Close();
951 delete fDisplaySocket;
953 fDisplaySocket = socket;
954 fDisplaySocket->SetOption(kNoBlock, 1);
955 TInetAddress adr = socket->GetInetAddress();
956 Info("CheckForConnections", "new display:\n %s (%s), port %d\n",
957 adr.GetHostName(), adr.GetHostAddress(), adr.GetPort());
961 for (Int_t iSocket = 0; iSocket < fSockets.GetEntriesFast(); iSocket++) {
962 socket = (TSocket*) fSockets[iSocket];
963 if (!socket) continue;
964 // remove finished client
966 if (socket->Recv(str, 255)) {
967 TString socketMessage(str);
968 if(socketMessage.CompareTo("Finished") == 0) {
969 TInetAddress adr = socket->GetInetAddress();
970 Info("CheckForConnections",
971 "disconnect finished client:\n %s (%s), port %d\n",
972 adr.GetHostName(), adr.GetHostAddress(), adr.GetPort());
973 delete fSockets.RemoveAt(iSocket);
977 if (!socket->IsValid()) {
978 // remove invalid sockets from the list
979 TInetAddress adr = socket->GetInetAddress();
980 Info("BroadcastHistos", "disconnect client:\n %s (%s), port %d\n",
981 adr.GetHostName(), adr.GetHostAddress(), adr.GetPort());
982 delete fSockets.RemoveAt(iSocket);
988 //_____________________________________________________________________________
989 void AliMonitorProcess::BroadcastHistos()
991 // send the monitor histograms to the clients
993 SetStatus(kBroadcasting);
994 TMessage message(kMESS_OBJECT);
995 message.WriteObject(fTopFolder);
997 for (Int_t iSocket = 0; iSocket < fSockets.GetEntriesFast(); iSocket++) {
998 TSocket* socket = (TSocket*) fSockets[iSocket];
999 if (!socket) continue;
1000 socket->SetOption(kNoBlock, 0);
1001 if (!socket->IsValid() || (socket->Send(message) < 0)) {
1002 // remove the socket from the list if there was an error
1003 TInetAddress adr = socket->GetInetAddress();
1004 Info("BroadcastHistos", "disconnect client:\n %s (%s), port %d\n",
1005 adr.GetHostName(), adr.GetHostAddress(), adr.GetPort());
1006 delete fSockets.RemoveAt(iSocket);
1008 socket->SetOption(kNoBlock, 1);
1011 fSockets.Compress();