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 "AliITSclustererV2.h"
39 #include "AliITSgeom.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"
55 #include <AliL3StandardIncludes.h>
56 #include <AliL3MemHandler.h>
57 #include <AliL3ClusterFitter.h>
58 #include <AliL3Fitter.h>
59 #include <AliL3Hough.h>
60 #include <AliL3HoughBaseTransformer.h>
61 #include <AliL3StandardIncludes.h>
62 #include <AliL3Track.h>
63 #include <AliL3TrackArray.h>
64 #include <AliL3Transform.h>
65 #include <AliL3Vertex.h>
66 #include <AliLevel3.h>
68 ClassImp(AliMonitorProcess)
71 const Int_t AliMonitorProcess::fgkPort = 9327;
74 //_____________________________________________________________________________
75 AliMonitorProcess::AliMonitorProcess(
76 #if ROOT_VERSION_CODE <= 199169 // 3.10/01
77 const char* /*alienHost*/,
79 const char* alienHost,
82 const char* selection,
83 const char* fileNameGalice):
84 fSelection(selection),
99 fWriteHistoList(kFALSE),
108 fDisplaySocket(NULL),
113 fInterruptHandler(NULL)
115 // initialize the monitoring process and the monitor histograms
117 fSelection = selection;
119 #if ROOT_VERSION_CODE <= 199169 // 3.10/01
120 fGrid = TGrid::Connect("alien", gSystem->Getenv("USER"));
122 fGrid = TGrid::Connect(alienHost, gSystem->Getenv("USER"));
124 if (!fGrid || fGrid->IsZombie() || !fGrid->IsConnected()) {
126 AliFatal("could not connect to alien");
128 #if ROOT_VERSION_CODE <= 199169 // 3.10/01
132 fRunLoader = AliRunLoader::Open(fileNameGalice);
133 if (!fRunLoader) AliFatal(Form("could not get run loader from file %s",
136 fRunLoader->CdGAFile();
137 fTPCParam = AliTPC::LoadTPCParam(gFile);
138 if (!fTPCParam) AliFatal("could not load TPC parameters");
140 fRunLoader->LoadgAlice();
141 gAlice = fRunLoader->GetAliRun();
142 if (!gAlice) AliFatal("no gAlice object found");
143 fITSgeom = (AliITSgeom*)gDirectory->Get("AliITSgeom");
144 if (!fITSgeom) AliFatal("could not load ITS geometry");
146 // Init TPC parameters for HLT
147 Bool_t isinit=AliL3Transform::Init(const_cast<char*>(fileNameGalice),kTRUE);
149 AliFatal("Could not create transform settings, please check log for error messages!");
152 fTopFolder = new TFolder("Monitor", "monitor histograms");
153 fTopFolder->SetOwner(kTRUE);
155 if (IsSelected("TPC")) fMonitors.Add(new AliMonitorTPC(fTPCParam));
156 if (IsSelected("ITS")) fMonitors.Add(new AliMonitorITS(fITSgeom));
157 if (IsSelected("V0s")) fMonitors.Add(new AliMonitorV0s);
158 if (IsSelected("HLTConfMap")) fMonitors.Add(new AliMonitorHLT(fTPCParam));
159 if (IsSelected("HLTHough")) fMonitors.Add(new AliMonitorHLTHough(fTPCParam));
161 for (Int_t iMonitor = 0; iMonitor < fMonitors.GetEntriesFast(); iMonitor++) {
162 ((AliMonitor*) fMonitors[iMonitor])->CreateHistos(fTopFolder);
165 TIterator* iFolder = fTopFolder->GetListOfFolders()->MakeIterator();
166 while (TFolder* folder = (TFolder*) iFolder->Next()) folder->SetOwner(kTRUE);
169 fFile = TFile::Open("monitor_tree.root", "RECREATE");
170 if (!fFile || !fFile->IsOpen()) {
171 AliFatal("could not open file for tree");
173 fTree = new TTree("MonitorTree", "tree for monitoring");
174 for (Int_t iMonitor = 0; iMonitor < fMonitors.GetEntriesFast(); iMonitor++) {
175 ((AliMonitor*) fMonitors[iMonitor])->CreateBranches(fTree);
179 fServerSocket = new TServerSocket(fgkPort, kTRUE);
180 fServerSocket->SetOption(kNoBlock, 1);
181 CheckForConnections();
183 fInterruptHandler = new AliMonitorInterruptHandler(this);
184 gSystem->AddSignalHandler(fInterruptHandler);
187 //_____________________________________________________________________________
188 AliMonitorProcess::AliMonitorProcess(const AliMonitorProcess& process) :
197 fLogicalFileName(""),
206 fWriteHistoList(kFALSE),
215 fDisplaySocket(NULL),
220 fInterruptHandler(NULL)
223 AliFatal("copy constructor not implemented");
226 //_____________________________________________________________________________
227 AliMonitorProcess& AliMonitorProcess::operator = (const AliMonitorProcess&
230 AliFatal("assignment operator not implemented");
234 //_____________________________________________________________________________
235 AliMonitorProcess::~AliMonitorProcess()
243 delete fServerSocket;
245 delete fDisplaySocket;
247 #if ROOT_VERSION_CODE <= 199169 // 3.10/01
254 gSystem->Unlink("monitor_tree.root");
259 gSystem->RemoveSignalHandler(fInterruptHandler);
260 delete fInterruptHandler;
264 //_____________________________________________________________________________
265 const char* AliMonitorProcess::GetRevision()
271 //_____________________________________________________________________________
272 void AliMonitorProcess::SetStatus(EStatus status)
274 // set the current status and process system events
277 gSystem->ProcessEvents();
281 //_____________________________________________________________________________
282 void AliMonitorProcess::Run()
284 // run the monitor process:
285 // check for a raw data file, process the raw data file and delete it
291 while (!CheckForNewFile()) {
292 CheckForConnections();
294 if (fStopping) break;
297 if (fStopping) break;
309 //_____________________________________________________________________________
310 void AliMonitorProcess::Stop()
312 // set the fStopping flag to terminate the monitor process after the current
313 // event was processed
315 if (GetStatus() != kStopped) fStopping = kTRUE;
319 //_____________________________________________________________________________
320 void AliMonitorProcess::ProcessFile(const char* fileName)
322 // create a file with monitor histograms for a single file
324 if (GetStatus() != kStopped) {
325 AliError("ProcessFile can not be called"
326 " while the monitor process is running");
330 fFileName = fileName;
331 Int_t nEventMin = fNEventsMin;
335 fNEventsMin = nEventMin;
340 //_____________________________________________________________________________
341 Bool_t AliMonitorProcess::CheckForNewFile()
343 // check whether a new file was registered in alien
345 #if ROOT_VERSION_CODE < ROOT_VERSION(5,0,0)
346 #if ROOT_VERSION_CODE <= 199169 // 3.10/01
347 TGridResult* result = fGrid->Ls();
351 sprintf(dirName, "%s/adc-%d", fAlienDir.Data(), datime.GetDate());
353 sprintf(findName, "*.root");
354 Grid_ResultHandle_t handle = fGrid->Find(dirName, findName);
356 AliError(Form("could not open alien directory %s",
360 TGridResult* result = fGrid->CreateGridResult(handle);
366 #if ROOT_VERSION_CODE <= 199169 // 3.10/01
367 while (const char* entry = result->Next()) {
369 while (Grid_Result_t* resultEntry = result->Next()) {
370 const char* entry = resultEntry->name.c_str();
372 if (strrchr(entry, '/')) entry = strrchr(entry, '/')+1;
373 // entry = host_date_time.root
374 TString entryCopy(entry);
375 char* p = const_cast<char*>(entryCopy.Data());
376 if (!strtok(p, "_") || !p) continue; // host name
377 char* dateStr = strtok(NULL, "_");
378 if (!dateStr || !p) continue;
379 char* timeStr = strtok(NULL, ".");
380 if (!timeStr || !p) continue;
381 Long_t date = atoi(dateStr);
382 Long_t time = atoi(timeStr);
384 if ((date > maxDate) || ((date == maxDate) && (time > maxTime))) {
392 if (maxDate < 0) return kFALSE; // no files found
393 if (fLogicalFileName.CompareTo(fileName) == 0) return kFALSE; // no new file
395 fLogicalFileName = fileName;
396 #if ROOT_VERSION_CODE <= 199169 // 3.10/01
397 result = fGrid->GetPhysicalFileNames(fLogicalFileName.Data());
398 fFileName = result->Next();
400 fileName = dirName + ("/" + fLogicalFileName);
401 handle = fGrid->GetPhysicalFileNames(fileName.Data());
403 AliError(Form("could not get physical file names for %s",
407 result = fGrid->CreateGridResult(handle);
409 Grid_Result_t* resultEntry = result->Next();
411 AliError(Form("could not get physical file names for %s",
415 fFileName = resultEntry->name2.c_str();
416 fFileName.ReplaceAll("castor:/", "rfio:/");
420 Error("CheckForNewFile", "needs to be ported to new TGrid");
425 //_____________________________________________________________________________
426 Bool_t AliMonitorProcess::ProcessFile()
428 // loop over all events in the raw data file, run the reconstruction
429 // and fill the monitor histograms
431 Int_t nEvents = GetNumberOfEvents(fFileName);
432 if (nEvents <= 0) return kFALSE;
433 AliDebug(1, Form("found %d event(s) in file %s",
434 nEvents, fFileName.Data()));
435 if (IsSelected("HLTConfMap")) CreateHLT(fFileName);
436 if (IsSelected("HLTHough")) CreateHLTHough(fFileName);
438 // loop over the events
439 for (Int_t iEvent = 0; iEvent < nEvents; iEvent++) {
440 CheckForConnections();
442 fRunLoader->SetEventNumber(0);
443 AliRawReaderRoot rawReader(fFileName, iEvent);
444 if (fStopping) break;
445 if (rawReader.GetRunNumber() != fRunNumber) {
448 fRunNumber = rawReader.GetRunNumber();
449 fEventNumber[0] = rawReader.GetEventId()[0];
450 fEventNumber[1] = rawReader.GetEventId()[1];
452 if (fStopping) break;
455 // monitor only central physics events
456 if (rawReader.GetType() != 7) continue;
457 if ((rawReader.GetAttributes()[0] & 0x02) == 0) continue;
458 AliInfo(Form("run: %d event: %d %d\n", rawReader.GetRunNumber(),
459 rawReader.GetEventId()[0], rawReader.GetEventId()[1]));
462 if (IsSelected("TPC")) {
463 CheckForConnections();
464 if (!ReconstructTPC(&rawReader, &esd)) return kFALSE;
465 if (fStopping) break;
467 if (IsSelected("ITS")) {
468 CheckForConnections();
469 if (!ReconstructITS(&rawReader, &esd)) return kFALSE;
470 if (fStopping) break;
472 if (IsSelected("V0s")) {
473 CheckForConnections();
474 if (!ReconstructV0s(&esd)) return kFALSE;
475 if (fStopping) break;
477 if (IsSelected("HLTConfMap")) {
478 CheckForConnections();
479 if (!ReconstructHLT(iEvent)) return kFALSE;
480 if (fStopping) break;
482 if (IsSelected("HLTHough")) {
483 CheckForConnections();
484 if (!ReconstructHLTHough(iEvent)) return kFALSE;
485 if (fStopping) break;
488 if (fDisplaySocket) fDisplaySocket->Send("new event");
490 AliDebug(1, "filling histograms...");
491 for (Int_t iMonitor = 0; iMonitor < fMonitors.GetEntriesFast(); iMonitor++) {
492 CheckForConnections();
494 ((AliMonitor*) fMonitors[iMonitor])->FillHistos(fRunLoader, &rawReader,
496 if (fStopping) break;
498 if (fStopping) break;
500 AliDebug(1, "updating histograms...");
501 CheckForConnections();
502 SetStatus(kUpdating);
503 TIterator* iFolder = fTopFolder->GetListOfFolders()->MakeIterator();
504 while (TFolder* folder = (TFolder*) iFolder->Next()) {
505 TIterator* iHisto = folder->GetListOfFolders()->MakeIterator();
506 while (AliMonitorPlot* histo = (AliMonitorPlot*) iHisto->Next()) {
512 if (fStopping) break;
514 AliDebug(1, "filling the tree...");
517 AliDebug(1, "broadcasting histograms...");
518 CheckForConnections();
522 if (fStopping) break;
537 //_____________________________________________________________________________
538 void AliMonitorProcess::Reset()
540 // write the current histograms to a file and reset them
542 if (fSubRunNumber == 0) fSubRunNumber++;
549 //_____________________________________________________________________________
550 UInt_t AliMonitorProcess::GetEventPeriodNumber() const
552 // get the period number from the event id
554 return (fEventNumber[1] >> 4);
557 //_____________________________________________________________________________
558 UInt_t AliMonitorProcess::GetEventOrbitNumber() const
560 // get the orbit number from the event id
562 return ((fEventNumber[1] & 0x000F) << 20) + (fEventNumber[0] >> 12);
565 //_____________________________________________________________________________
566 UInt_t AliMonitorProcess::GetEventBunchNumber() const
568 // get the bunch number from the event id
570 return (fEventNumber[0] % 0x0FFF);
573 //_____________________________________________________________________________
574 Int_t AliMonitorProcess::GetNumberOfEvents(const char* fileName) const
576 // determine the number of events in the given raw data file
580 TFile* file = TFile::Open(fileName);
581 if (!file || !file->IsOpen()) {
582 AliError(Form("could not open file %s", fileName));
583 if (file) delete file;
587 TTree* tree = (TTree*) file->Get("RAW");
589 AliError("could not find tree with raw data");
591 nEvents = (Int_t) tree->GetEntries();
599 //_____________________________________________________________________________
600 Bool_t AliMonitorProcess::ReconstructTPC(AliRawReader* rawReader, AliESD* esd)
602 // find TPC clusters and tracks
606 AliLoader* tpcLoader = fRunLoader->GetLoader("TPCLoader");
608 AliError("no TPC loader found");
611 gSystem->Unlink("TPC.RecPoints.root");
614 AliDebug(1, "reconstructing clusters...");
615 tpcLoader->LoadRecPoints("recreate");
616 AliTPCclustererMI clusterer(fTPCParam);
617 tpcLoader->MakeRecPointsContainer();
618 clusterer.SetOutput(tpcLoader->TreeR());
619 clusterer.Digits2Clusters(rawReader);
620 tpcLoader->WriteRecPoints("OVERWRITE");
623 AliDebug(1, "reconstructing tracks...");
624 AliTPCtrackerMI tracker(fTPCParam);
625 tracker.LoadClusters(tpcLoader->TreeR());
626 tracker.Clusters2Tracks(esd);
627 tracker.UnloadClusters();
628 tpcLoader->UnloadRecPoints();
633 //_____________________________________________________________________________
634 Bool_t AliMonitorProcess::ReconstructITS(AliRawReader* rawReader, AliESD* esd)
636 // find ITS clusters and tracks
640 AliLoader* itsLoader = fRunLoader->GetLoader("ITSLoader");
642 AliError("no ITS loader found");
645 gSystem->Unlink("ITS.RecPoints.root");
648 AliDebug(1, "reconstructing clusters...");
649 itsLoader->LoadRecPoints("recreate");
650 AliITSclustererV2 clusterer(fITSgeom);
651 itsLoader->MakeRecPointsContainer();
652 clusterer.Digits2Clusters(rawReader);
655 AliDebug(1, "reconstructing tracks...");
656 AliITStrackerV2 tracker(fITSgeom);
657 tracker.LoadClusters(itsLoader->TreeR());
658 tracker.Clusters2Tracks(esd);
659 tracker.UnloadClusters();
661 itsLoader->UnloadRecPoints();
665 //_____________________________________________________________________________
666 Bool_t AliMonitorProcess::ReconstructV0s(AliESD* esd)
673 AliDebug(1, "reconstructing V0s...");
674 AliV0vertexer vertexer;
676 esd->GetVertex()->GetXYZ(vtx);
677 vertexer.SetVertex(vtx);
678 vertexer.Tracks2V0vertices(esd);
683 //_____________________________________________________________________________
684 void AliMonitorProcess::CreateHLT(const char* fileName)
687 // create the HLT (Level3) object
689 if (fHLT) delete fHLT;
692 strcpy(name, fileName);
693 fHLT = new AliLevel3(name);
694 fHLT->Init("./", AliLevel3::kRaw, 1);
696 fHLT->SetClusterFinderParam(-1, -1, kTRUE);
698 Int_t phiSegments = 50;
699 Int_t etaSegments = 100;
700 Int_t trackletlength = 3;
701 Int_t tracklength = 20;//40 or 5
702 Int_t rowscopetracklet = 2;
703 Int_t rowscopetrack = 10;
704 Double_t minPtFit = 0;
705 Double_t maxangle = 0.1745;
706 Double_t goodDist = 5;
707 Double_t maxphi = 0.1;
708 Double_t maxeta = 0.1;
709 Double_t hitChi2Cut = 15;//100 or 15
710 Double_t goodHitChi2 = 5;//20 or 5
711 Double_t trackChi2Cut = 10;//50 or 10
712 fHLT->SetTrackerParam(phiSegments, etaSegments,
713 trackletlength, tracklength,
714 rowscopetracklet, rowscopetrack,
715 minPtFit, maxangle, goodDist, hitChi2Cut,
716 goodHitChi2, trackChi2Cut, 50, maxphi, maxeta, kTRUE);
718 fHLT->WriteFiles("./hlt/");
721 //_____________________________________________________________________________
722 void AliMonitorProcess::CreateHLTHough(const char* fileName)
725 // create the HLT Hough transform (L3Hough) object
727 if (fHLTHough) delete fHLTHough;
730 strcpy(name, fileName);
732 fHLTHough = new AliL3Hough();
733 fHLTHough->SetThreshold(4);
734 fHLTHough->SetTransformerParams(140,150,0.5,-1);
735 fHLTHough->SetPeakThreshold(9000,-1);// or 6000
736 fHLTHough->Init("./", kFALSE, 50, kFALSE,0,name);
737 fHLTHough->SetAddHistograms();
738 // fHLTHough->GetMaxFinder()->SetThreshold(14000);
742 //_____________________________________________________________________________
743 Bool_t AliMonitorProcess::ReconstructHLT(Int_t iEvent)
745 // run the HLT cluster and track finder
749 gSystem->Exec("rm -rf hlt");
750 gSystem->MakeDirectory("hlt");
751 if (!fHLT) return kFALSE;
753 fHLT->ProcessEvent(0, 35, iEvent);
755 // remove the event number from the file names
757 sprintf(command, "rename points_%d points hlt/*.raw", iEvent);
758 gSystem->Exec(command);
759 sprintf(command, "rename tracks_tr_%d tracks_tr hlt/*.raw", iEvent);
760 gSystem->Exec(command);
761 sprintf(command, "rename tracks_gl_%d tracks_gl hlt/*.raw", iEvent);
762 gSystem->Exec(command);
763 sprintf(command, "rename tracks_%d tracks hlt/*.raw", iEvent);
764 gSystem->Exec(command);
768 //_____________________________________________________________________________
769 Bool_t AliMonitorProcess::ReconstructHLTHough(Int_t iEvent)
771 // run the HLT Hough transformer
775 gSystem->Exec("rm -rf hlt/hough");
776 gSystem->MakeDirectory("hlt/hough");
777 gSystem->Exec("rm -rf hlt/fitter");
778 gSystem->MakeDirectory("hlt/fitter");
779 if (!fHLTHough) return kFALSE;
781 // fHLTHough->Process(0, 35);
782 // Loop over TPC sectors and process the data
783 for(Int_t i=0; i<=35; i++)
785 fHLTHough->ReadData(i,iEvent);
786 fHLTHough->Transform();
787 // if(fHLTHough->fAddHistograms)
788 fHLTHough->AddAllHistograms();
789 fHLTHough->FindTrackCandidates();
790 fHLTHough->AddTracks();
792 fHLTHough->WriteTracks("./hlt/hough");
794 // Run cluster fitter
795 AliL3ClusterFitter *fitter = new AliL3ClusterFitter("./hlt");
797 // Set debug flag for the cluster fitter
800 // Setting fitter parameters
801 fitter->SetInnerWidthFactor(1,1.5);
802 fitter->SetOuterWidthFactor(1,1.5);
803 fitter->SetNmaxOverlaps(5);
805 //fitter->SetChiSqMax(5,kFALSE); //isolated clusters
806 fitter->SetChiSqMax(5,kTRUE); //overlapping clusters
808 Int_t rowrange[2] = {0,AliL3Transform::GetNRows()-1};
810 // Takes input from global hough tracks produced by HT
811 fitter->LoadSeeds(rowrange,kFALSE,iEvent);
815 for(Int_t islice = 0; islice <= 35; islice++)
817 for(Int_t ipatch = 0; ipatch < AliL3Transform::GetNPatches(); ipatch++)
820 fHLTHough->GetMemHandler(ipatch)->Free();
821 fHLTHough->GetMemHandler(ipatch)->Init(islice,ipatch);
822 AliL3DigitRowData *digits = (AliL3DigitRowData *)fHLTHough->GetMemHandler(ipatch)->AliAltroDigits2Memory(ndigits,iEvent);
824 fitter->Init(islice,ipatch);
825 fitter->SetInputData(digits);
826 fitter->FindClusters();
827 fitter->WriteClusters();
831 // Refit of the clusters
833 //The seeds are the input tracks from circle HT
834 AliL3TrackArray *tracks = fitter->GetSeeds();
835 AliL3Fitter *ft = new AliL3Fitter(&vertex,1);
837 ft->LoadClusters("./hlt/fitter/",iEvent,kFALSE);
838 for(Int_t i=0; i<tracks->GetNTracks(); i++)
840 AliL3Track *track = tracks->GetCheckedTrack(i);
842 if(track->GetNHits() < 20) continue;
843 ft->SortTrackClusters(track);
845 track->UpdateToFirstPoint();
849 //Write the final tracks
850 fitter->WriteTracks(20);
854 // remove the event number from the file names
856 sprintf(command, "rename tracks_%d tracks hlt/hough/*.raw", iEvent);
857 gSystem->Exec(command);
858 sprintf(command, "rename tracks_%d tracks hlt/fitter/*.raw", iEvent);
859 gSystem->Exec(command);
860 sprintf(command, "rename points_%d points hlt/fitter/*.raw", iEvent);
861 gSystem->Exec(command);
865 //_____________________________________________________________________________
866 Bool_t AliMonitorProcess::WriteHistos()
868 // write the monitor tree and the monitor histograms to the file
869 // "monitor_<run number>[_<sub_run_number>].root"
870 // if at least fNEventsMin events were monitored
874 // rename tree file and create a new one
881 sprintf(fileName, "monitor_tree_%d.root", fRunNumber);
882 if (fSubRunNumber > 0) {
883 sprintf(fileName, "monitor_tree_%d_%d.root", fRunNumber, fSubRunNumber);
885 if (fNEvents < fNEventsMin) {
886 gSystem->Unlink("monitor_tree.root");
888 gSystem->Rename("monitor_tree.root", fileName);
891 fFile = TFile::Open("monitor_tree.root", "RECREATE");
892 if (!fFile || !fFile->IsOpen()) {
893 AliFatal("could not open file for tree");
895 fTree = new TTree("MonitorTree", "tree for monitoring");
896 for (Int_t iMonitor = 0; iMonitor < fMonitors.GetEntriesFast(); iMonitor++) {
897 ((AliMonitor*) fMonitors[iMonitor])->CreateBranches(fTree);
901 // write the histograms
902 if (fNEvents < fNEventsMin) return kTRUE;
904 if (!fWriteHistoList) {
905 TIterator* iFolder = fTopFolder->GetListOfFolders()->MakeIterator();
906 while (TFolder* folder = (TFolder*) iFolder->Next()) {
907 TIterator* iHisto = folder->GetListOfFolders()->MakeIterator();
908 while (AliMonitorPlot* histo = (AliMonitorPlot*) iHisto->Next()) {
916 Bool_t result = kTRUE;
917 sprintf(fileName, "monitor_%d.root", fRunNumber);
918 if (fSubRunNumber > 0) {
919 sprintf(fileName, "monitor_%d_%d.root", fRunNumber, fSubRunNumber);
921 TFile* file = TFile::Open(fileName, "recreate");
922 if (!file || !file->IsOpen()) {
923 AliError(Form("could not open file %s", fileName));
929 if (file) delete file;
934 //_____________________________________________________________________________
935 void AliMonitorProcess::StartNewRun()
937 // reset the histograms for a new run
939 SetStatus(kResetting);
940 TIterator* iFolder = fTopFolder->GetListOfFolders()->MakeIterator();
941 while (TFolder* folder = (TFolder*) iFolder->Next()) {
942 TIterator* iHisto = folder->GetListOfFolders()->MakeIterator();
943 while (AliMonitorPlot* histo = (AliMonitorPlot*) iHisto->Next()) {
954 //_____________________________________________________________________________
955 void AliMonitorProcess::CheckForConnections()
957 // check if new clients want to connect and add them to the list of sockets
960 while ((socket = fServerSocket->Accept()) != (TSocket*)-1) {
961 socket->SetOption(kNoBlock, 1);
962 char socketType[256];
963 if (socket->Recv(socketType, 255) <= 0) {
964 gSystem->Sleep(1000);
965 if (socket->Recv(socketType, 255) <= 0) {
966 TInetAddress adr = socket->GetInetAddress();
967 AliError(Form("no socket type received - "
968 "disconnect client: %s (%s), port %d",
969 adr.GetHostName(), adr.GetHostAddress(), adr.GetPort()));
974 if (strcmp(socketType, "client") == 0) {
975 fSockets.Add(socket);
976 TInetAddress adr = socket->GetInetAddress();
977 AliInfo(Form("new client: %s (%s), port %d",
978 adr.GetHostName(), adr.GetHostAddress(), adr.GetPort()));
979 if (fNEvents > 0) BroadcastHistos(socket);
980 } else if (strcmp(socketType, "display") == 0) {
981 if (fDisplaySocket) {
982 fDisplaySocket->Close();
983 delete fDisplaySocket;
985 fDisplaySocket = socket;
986 TInetAddress adr = socket->GetInetAddress();
987 AliInfo(Form("new display: %s (%s), port %d",
988 adr.GetHostName(), adr.GetHostAddress(), adr.GetPort()));
990 TInetAddress adr = socket->GetInetAddress();
991 AliError(Form("unknown socket type %s - "
992 "disconnect client: %s (%s), port %d", socketType,
993 adr.GetHostName(), adr.GetHostAddress(), adr.GetPort()));
999 // remove finished or invalid clients
1000 for (Int_t iSocket = 0; iSocket < fSockets.GetEntriesFast(); iSocket++) {
1001 socket = (TSocket*) fSockets[iSocket];
1002 if (!socket) continue;
1003 char controlMessage[256];
1004 if (socket->Recv(controlMessage, 255)) {
1005 if (strcmp(controlMessage, "disconnect") == 0) {
1006 TInetAddress adr = socket->GetInetAddress();
1007 AliInfo(Form("disconnect client: %s (%s), port %d",
1008 adr.GetHostName(), adr.GetHostAddress(), adr.GetPort()));
1009 delete fSockets.RemoveAt(iSocket);
1013 if (!socket->IsValid()) {
1014 // remove invalid sockets from the list
1015 TInetAddress adr = socket->GetInetAddress();
1016 AliError(Form("disconnect invalid client: %s (%s), port %d",
1017 adr.GetHostName(), adr.GetHostAddress(), adr.GetPort()));
1018 delete fSockets.RemoveAt(iSocket);
1021 fSockets.Compress();
1024 //_____________________________________________________________________________
1025 void AliMonitorProcess::BroadcastHistos(TSocket* toSocket)
1027 // send the monitor histograms to the clients
1029 SetStatus(kBroadcasting);
1030 TMessage message(kMESS_OBJECT);
1031 message.WriteObject(fTopFolder);
1033 for (Int_t iSocket = 0; iSocket < fSockets.GetEntriesFast(); iSocket++) {
1034 TSocket* socket = (TSocket*) fSockets[iSocket];
1035 if (!socket) continue;
1036 if (toSocket && (socket != toSocket)) continue;
1038 // send control message
1039 if (!socket->IsValid() || (socket->Send("histograms") <= 0)) {
1040 TInetAddress adr = socket->GetInetAddress();
1041 AliError(Form("connection to client failed - "
1042 "disconnect client: %s (%s), port %d",
1043 adr.GetHostName(), adr.GetHostAddress(), adr.GetPort()));
1044 delete fSockets.RemoveAt(iSocket);
1047 // receive control message
1048 char controlMessage[256];
1049 Int_t result = socket->Recv(controlMessage, 255);
1051 gSystem->Sleep(1000); // wait one second and try again
1052 result = socket->Recv(controlMessage, 255);
1055 TInetAddress adr = socket->GetInetAddress();
1056 AliError(Form("no response from client - "
1057 "disconnect client: %s (%s), port %d",
1058 adr.GetHostName(), adr.GetHostAddress(), adr.GetPort()));
1059 delete fSockets.RemoveAt(iSocket);
1062 if (strcmp(controlMessage, "ok") != 0) {
1063 TInetAddress adr = socket->GetInetAddress();
1064 AliError(Form("no \"ok\" message from client - "
1065 "disconnect client: %s (%s), port %d",
1066 adr.GetHostName(), adr.GetHostAddress(), adr.GetPort()));
1067 delete fSockets.RemoveAt(iSocket);
1071 socket->SetOption(kNoBlock, 0);
1072 if (socket->Send(message) < 0) {
1073 // remove the socket from the list if there was an error
1074 TInetAddress adr = socket->GetInetAddress();
1075 AliError(Form("sending histograms failed - "
1076 "disconnect client: %s (%s), port %d",
1077 adr.GetHostName(), adr.GetHostAddress(), adr.GetPort()));
1078 delete fSockets.RemoveAt(iSocket);
1080 gSystem->Sleep(100);
1081 socket->SetOption(kNoBlock, 1);
1084 fSockets.Compress();
1088 //_____________________________________________________________________________
1089 AliMonitorProcess::AliMonitorInterruptHandler::AliMonitorInterruptHandler
1090 (AliMonitorProcess* process):
1091 TSignalHandler(kSigUser1, kFALSE),
1094 // constructor: set process
1097 //_____________________________________________________________________________
1098 AliMonitorProcess::AliMonitorInterruptHandler::AliMonitorInterruptHandler
1099 (const AliMonitorInterruptHandler& handler):
1100 TSignalHandler(handler)
1104 AliFatal("copy constructor not implemented");
1107 //_____________________________________________________________________________
1108 AliMonitorProcess::AliMonitorInterruptHandler&
1109 AliMonitorProcess::AliMonitorInterruptHandler::operator =
1110 (const AliMonitorInterruptHandler& /*handler*/)
1112 // assignment operator
1114 AliFatal("assignment operator not implemented");
1118 //_____________________________________________________________________________
1119 Bool_t AliMonitorProcess::AliMonitorInterruptHandler::Notify()
1121 // interrupt signal -> stop process
1123 AliInfo("the monitoring process will be stopped.");