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>
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"
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 AliITS* its = (AliITS*) gAlice->GetModule("ITS");
144 if (!its) AliFatal("no ITS detector found");
145 fITSgeom = its->GetITSgeom();
146 if (!fITSgeom) AliFatal("could not load ITS geometry");
148 // Init TPC parameters for HLT
149 Bool_t isinit=AliL3Transform::Init(const_cast<char*>(fileNameGalice),kTRUE);
151 AliFatal("Could not create transform settings, please check log for error messages!");
154 fTopFolder = new TFolder("Monitor", "monitor histograms");
155 fTopFolder->SetOwner(kTRUE);
157 if (IsSelected("TPC")) fMonitors.Add(new AliMonitorTPC(fTPCParam));
158 if (IsSelected("ITS")) fMonitors.Add(new AliMonitorITS(fITSgeom));
159 if (IsSelected("V0s")) fMonitors.Add(new AliMonitorV0s);
160 if (IsSelected("HLTConfMap")) fMonitors.Add(new AliMonitorHLT(fTPCParam));
161 if (IsSelected("HLTHough")) fMonitors.Add(new AliMonitorHLTHough(fTPCParam));
163 for (Int_t iMonitor = 0; iMonitor < fMonitors.GetEntriesFast(); iMonitor++) {
164 ((AliMonitor*) fMonitors[iMonitor])->CreateHistos(fTopFolder);
167 TIterator* iFolder = fTopFolder->GetListOfFolders()->MakeIterator();
168 while (TFolder* folder = (TFolder*) iFolder->Next()) folder->SetOwner(kTRUE);
171 fFile = TFile::Open("monitor_tree.root", "RECREATE");
172 if (!fFile || !fFile->IsOpen()) {
173 AliFatal("could not open file for tree");
175 fTree = new TTree("MonitorTree", "tree for monitoring");
176 for (Int_t iMonitor = 0; iMonitor < fMonitors.GetEntriesFast(); iMonitor++) {
177 ((AliMonitor*) fMonitors[iMonitor])->CreateBranches(fTree);
181 fServerSocket = new TServerSocket(fgkPort, kTRUE);
182 fServerSocket->SetOption(kNoBlock, 1);
183 CheckForConnections();
185 fInterruptHandler = new AliMonitorInterruptHandler(this);
186 gSystem->AddSignalHandler(fInterruptHandler);
189 //_____________________________________________________________________________
190 AliMonitorProcess::AliMonitorProcess(const AliMonitorProcess& process) :
199 fLogicalFileName(""),
208 fWriteHistoList(kFALSE),
217 fDisplaySocket(NULL),
222 fInterruptHandler(NULL)
225 AliFatal("copy constructor not implemented");
228 //_____________________________________________________________________________
229 AliMonitorProcess& AliMonitorProcess::operator = (const AliMonitorProcess&
232 AliFatal("assignment operator not implemented");
236 //_____________________________________________________________________________
237 AliMonitorProcess::~AliMonitorProcess()
245 delete fServerSocket;
247 delete fDisplaySocket;
249 #if ROOT_VERSION_CODE <= 199169 // 3.10/01
256 gSystem->Unlink("monitor_tree.root");
261 gSystem->RemoveSignalHandler(fInterruptHandler);
262 delete fInterruptHandler;
266 //_____________________________________________________________________________
267 const char* AliMonitorProcess::GetRevision()
273 //_____________________________________________________________________________
274 void AliMonitorProcess::SetStatus(EStatus status)
276 // set the current status and process system events
279 gSystem->ProcessEvents();
283 //_____________________________________________________________________________
284 void AliMonitorProcess::Run()
286 // run the monitor process:
287 // check for a raw data file, process the raw data file and delete it
293 while (!CheckForNewFile()) {
294 CheckForConnections();
296 if (fStopping) break;
299 if (fStopping) break;
311 //_____________________________________________________________________________
312 void AliMonitorProcess::Stop()
314 // set the fStopping flag to terminate the monitor process after the current
315 // event was processed
317 if (GetStatus() != kStopped) fStopping = kTRUE;
321 //_____________________________________________________________________________
322 void AliMonitorProcess::ProcessFile(const char* fileName)
324 // create a file with monitor histograms for a single file
326 if (GetStatus() != kStopped) {
327 AliError("ProcessFile can not be called"
328 " while the monitor process is running");
332 fFileName = fileName;
333 Int_t nEventMin = fNEventsMin;
337 fNEventsMin = nEventMin;
342 //_____________________________________________________________________________
343 Bool_t AliMonitorProcess::CheckForNewFile()
345 // check whether a new file was registered in alien
347 #if ROOT_VERSION_CODE <= 199169 // 3.10/01
348 TGridResult* result = fGrid->Ls();
352 sprintf(dirName, "%s/adc-%d", fAlienDir.Data(), datime.GetDate());
354 sprintf(findName, "*.root");
355 Grid_ResultHandle_t handle = fGrid->Find(dirName, findName);
357 AliError(Form("could not open alien directory %s",
361 TGridResult* result = fGrid->CreateGridResult(handle);
367 #if ROOT_VERSION_CODE <= 199169 // 3.10/01
368 while (const char* entry = result->Next()) {
370 while (Grid_Result_t* resultEntry = result->Next()) {
371 const char* entry = resultEntry->name.c_str();
373 if (strrchr(entry, '/')) entry = strrchr(entry, '/')+1;
374 // entry = host_date_time.root
375 TString entryCopy(entry);
376 char* p = const_cast<char*>(entryCopy.Data());
377 if (!strtok(p, "_") || !p) continue; // host name
378 char* dateStr = strtok(NULL, "_");
379 if (!dateStr || !p) continue;
380 char* timeStr = strtok(NULL, ".");
381 if (!timeStr || !p) continue;
382 Long_t date = atoi(dateStr);
383 Long_t time = atoi(timeStr);
385 if ((date > maxDate) || ((date == maxDate) && (time > maxTime))) {
393 if (maxDate < 0) return kFALSE; // no files found
394 if (fLogicalFileName.CompareTo(fileName) == 0) return kFALSE; // no new file
396 fLogicalFileName = fileName;
397 #if ROOT_VERSION_CODE <= 199169 // 3.10/01
398 result = fGrid->GetPhysicalFileNames(fLogicalFileName.Data());
399 fFileName = result->Next();
401 fileName = dirName + ("/" + fLogicalFileName);
402 handle = fGrid->GetPhysicalFileNames(fileName.Data());
404 AliError(Form("could not get physical file names for %s",
408 result = fGrid->CreateGridResult(handle);
410 Grid_Result_t* resultEntry = result->Next();
412 AliError(Form("could not get physical file names for %s",
416 fFileName = resultEntry->name2.c_str();
417 fFileName.ReplaceAll("castor:/", "rfio:/");
424 //_____________________________________________________________________________
425 Bool_t AliMonitorProcess::ProcessFile()
427 // loop over all events in the raw data file, run the reconstruction
428 // and fill the monitor histograms
430 Int_t nEvents = GetNumberOfEvents(fFileName);
431 if (nEvents <= 0) return kFALSE;
432 AliDebug(1, Form("found %d event(s) in file %s",
433 nEvents, fFileName.Data()));
434 if (IsSelected("HLTConfMap")) CreateHLT(fFileName);
435 if (IsSelected("HLTHough")) CreateHLTHough(fFileName);
437 // loop over the events
438 for (Int_t iEvent = 0; iEvent < nEvents; iEvent++) {
439 CheckForConnections();
441 fRunLoader->SetEventNumber(0);
442 AliRawReaderRoot rawReader(fFileName, iEvent);
443 if (fStopping) break;
444 if (rawReader.GetRunNumber() != fRunNumber) {
447 fRunNumber = rawReader.GetRunNumber();
448 fEventNumber[0] = rawReader.GetEventId()[0];
449 fEventNumber[1] = rawReader.GetEventId()[1];
451 if (fStopping) break;
454 // monitor only central physics events
455 if (rawReader.GetType() != 7) continue;
456 if ((rawReader.GetAttributes()[0] & 0x02) == 0) continue;
457 AliInfo(Form("run: %d event: %d %d\n", rawReader.GetRunNumber(),
458 rawReader.GetEventId()[0], rawReader.GetEventId()[1]));
461 if (IsSelected("TPC")) {
462 CheckForConnections();
463 if (!ReconstructTPC(&rawReader, &esd)) return kFALSE;
464 if (fStopping) break;
466 if (IsSelected("ITS")) {
467 CheckForConnections();
468 if (!ReconstructITS(&rawReader, &esd)) return kFALSE;
469 if (fStopping) break;
471 if (IsSelected("V0s")) {
472 CheckForConnections();
473 if (!ReconstructV0s(&esd)) return kFALSE;
474 if (fStopping) break;
476 if (IsSelected("HLTConfMap")) {
477 CheckForConnections();
478 if (!ReconstructHLT(iEvent)) return kFALSE;
479 if (fStopping) break;
481 if (IsSelected("HLTHough")) {
482 CheckForConnections();
483 if (!ReconstructHLTHough(iEvent)) return kFALSE;
484 if (fStopping) break;
487 if (fDisplaySocket) fDisplaySocket->Send("new event");
489 AliDebug(1, "filling histograms...");
490 for (Int_t iMonitor = 0; iMonitor < fMonitors.GetEntriesFast(); iMonitor++) {
491 CheckForConnections();
493 ((AliMonitor*) fMonitors[iMonitor])->FillHistos(fRunLoader, &rawReader,
495 if (fStopping) break;
497 if (fStopping) break;
499 AliDebug(1, "updating histograms...");
500 CheckForConnections();
501 SetStatus(kUpdating);
502 TIterator* iFolder = fTopFolder->GetListOfFolders()->MakeIterator();
503 while (TFolder* folder = (TFolder*) iFolder->Next()) {
504 TIterator* iHisto = folder->GetListOfFolders()->MakeIterator();
505 while (AliMonitorPlot* histo = (AliMonitorPlot*) iHisto->Next()) {
511 if (fStopping) break;
513 AliDebug(1, "filling the tree...");
516 AliDebug(1, "broadcasting histograms...");
517 CheckForConnections();
521 if (fStopping) break;
536 //_____________________________________________________________________________
537 void AliMonitorProcess::Reset()
539 // write the current histograms to a file and reset them
541 if (fSubRunNumber == 0) fSubRunNumber++;
548 //_____________________________________________________________________________
549 UInt_t AliMonitorProcess::GetEventPeriodNumber() const
551 // get the period number from the event id
553 return (fEventNumber[1] >> 4);
556 //_____________________________________________________________________________
557 UInt_t AliMonitorProcess::GetEventOrbitNumber() const
559 // get the orbit number from the event id
561 return ((fEventNumber[1] & 0x000F) << 20) + (fEventNumber[0] >> 12);
564 //_____________________________________________________________________________
565 UInt_t AliMonitorProcess::GetEventBunchNumber() const
567 // get the bunch number from the event id
569 return (fEventNumber[0] % 0x0FFF);
572 //_____________________________________________________________________________
573 Int_t AliMonitorProcess::GetNumberOfEvents(const char* fileName) const
575 // determine the number of events in the given raw data file
579 TFile* file = TFile::Open(fileName);
580 if (!file || !file->IsOpen()) {
581 AliError(Form("could not open file %s", fileName));
582 if (file) delete file;
586 TTree* tree = (TTree*) file->Get("RAW");
588 AliError("could not find tree with raw data");
590 nEvents = (Int_t) tree->GetEntries();
598 //_____________________________________________________________________________
599 Bool_t AliMonitorProcess::ReconstructTPC(AliRawReader* rawReader, AliESD* esd)
601 // find TPC clusters and tracks
605 AliLoader* tpcLoader = fRunLoader->GetLoader("TPCLoader");
607 AliError("no TPC loader found");
610 gSystem->Unlink("TPC.RecPoints.root");
613 AliDebug(1, "reconstructing clusters...");
614 tpcLoader->LoadRecPoints("recreate");
615 AliTPCclustererMI clusterer(fTPCParam);
616 tpcLoader->MakeRecPointsContainer();
617 clusterer.SetOutput(tpcLoader->TreeR());
618 clusterer.Digits2Clusters(rawReader);
619 tpcLoader->WriteRecPoints("OVERWRITE");
622 AliDebug(1, "reconstructing tracks...");
623 AliTPCtrackerMI tracker(fTPCParam);
624 tracker.LoadClusters(tpcLoader->TreeR());
625 tracker.Clusters2Tracks(esd);
626 tracker.UnloadClusters();
627 tpcLoader->UnloadRecPoints();
632 //_____________________________________________________________________________
633 Bool_t AliMonitorProcess::ReconstructITS(AliRawReader* rawReader, AliESD* esd)
635 // find ITS clusters and tracks
639 AliLoader* itsLoader = fRunLoader->GetLoader("ITSLoader");
641 AliError("no ITS loader found");
644 gSystem->Unlink("ITS.RecPoints.root");
647 AliDebug(1, "reconstructing clusters...");
648 itsLoader->LoadRecPoints("recreate");
649 AliITSclustererV2 clusterer(fITSgeom);
650 itsLoader->MakeRecPointsContainer();
651 clusterer.Digits2Clusters(rawReader);
654 AliDebug(1, "reconstructing tracks...");
655 AliITStrackerV2 tracker(fITSgeom);
656 tracker.LoadClusters(itsLoader->TreeR());
657 tracker.Clusters2Tracks(esd);
658 tracker.UnloadClusters();
660 itsLoader->UnloadRecPoints();
664 //_____________________________________________________________________________
665 Bool_t AliMonitorProcess::ReconstructV0s(AliESD* esd)
672 AliDebug(1, "reconstructing V0s...");
673 AliV0vertexer vertexer;
675 esd->GetVertex()->GetXYZ(vtx);
676 vertexer.SetVertex(vtx);
677 vertexer.Tracks2V0vertices(esd);
682 //_____________________________________________________________________________
683 void AliMonitorProcess::CreateHLT(const char* fileName)
686 // create the HLT (Level3) object
688 if (fHLT) delete fHLT;
691 strcpy(name, fileName);
692 fHLT = new AliLevel3(name);
693 fHLT->Init("./", AliLevel3::kRaw, 1);
695 fHLT->SetClusterFinderParam(-1, -1, kTRUE);
697 Int_t phiSegments = 50;
698 Int_t etaSegments = 100;
699 Int_t trackletlength = 3;
700 Int_t tracklength = 20;//40 or 5
701 Int_t rowscopetracklet = 2;
702 Int_t rowscopetrack = 10;
703 Double_t minPtFit = 0;
704 Double_t maxangle = 0.1745;
705 Double_t goodDist = 5;
706 Double_t maxphi = 0.1;
707 Double_t maxeta = 0.1;
708 Double_t hitChi2Cut = 15;//100 or 15
709 Double_t goodHitChi2 = 5;//20 or 5
710 Double_t trackChi2Cut = 10;//50 or 10
711 fHLT->SetTrackerParam(phiSegments, etaSegments,
712 trackletlength, tracklength,
713 rowscopetracklet, rowscopetrack,
714 minPtFit, maxangle, goodDist, hitChi2Cut,
715 goodHitChi2, trackChi2Cut, 50, maxphi, maxeta, kTRUE);
717 fHLT->WriteFiles("./hlt/");
720 //_____________________________________________________________________________
721 void AliMonitorProcess::CreateHLTHough(const char* fileName)
724 // create the HLT Hough transform (L3Hough) object
726 if (fHLTHough) delete fHLTHough;
729 strcpy(name, fileName);
731 fHLTHough = new AliL3Hough();
732 fHLTHough->SetThreshold(4);
733 fHLTHough->SetTransformerParams(140,150,0.5,-1);
734 fHLTHough->SetPeakThreshold(9000,-1);// or 6000
735 fHLTHough->Init("./", kFALSE, 50, kFALSE,0,name);
736 fHLTHough->SetAddHistograms();
737 // fHLTHough->GetMaxFinder()->SetThreshold(14000);
741 //_____________________________________________________________________________
742 Bool_t AliMonitorProcess::ReconstructHLT(Int_t iEvent)
744 // run the HLT cluster and track finder
748 gSystem->Exec("rm -rf hlt");
749 gSystem->MakeDirectory("hlt");
750 if (!fHLT) return kFALSE;
752 fHLT->ProcessEvent(0, 35, iEvent);
754 // remove the event number from the file names
756 sprintf(command, "rename points_%d points hlt/*.raw", iEvent);
757 gSystem->Exec(command);
758 sprintf(command, "rename tracks_tr_%d tracks_tr hlt/*.raw", iEvent);
759 gSystem->Exec(command);
760 sprintf(command, "rename tracks_gl_%d tracks_gl hlt/*.raw", iEvent);
761 gSystem->Exec(command);
762 sprintf(command, "rename tracks_%d tracks hlt/*.raw", iEvent);
763 gSystem->Exec(command);
767 //_____________________________________________________________________________
768 Bool_t AliMonitorProcess::ReconstructHLTHough(Int_t iEvent)
770 // run the HLT Hough transformer
774 gSystem->Exec("rm -rf hlt/hough");
775 gSystem->MakeDirectory("hlt/hough");
776 gSystem->Exec("rm -rf hlt/fitter");
777 gSystem->MakeDirectory("hlt/fitter");
778 if (!fHLTHough) return kFALSE;
780 // fHLTHough->Process(0, 35);
781 // Loop over TPC sectors and process the data
782 for(Int_t i=0; i<=35; i++)
784 fHLTHough->ReadData(i,iEvent);
785 fHLTHough->Transform();
786 // if(fHLTHough->fAddHistograms)
787 fHLTHough->AddAllHistograms();
788 fHLTHough->FindTrackCandidates();
789 fHLTHough->AddTracks();
791 fHLTHough->WriteTracks("./hlt/hough");
793 // Run cluster fitter
794 AliL3ClusterFitter *fitter = new AliL3ClusterFitter("./hlt");
796 // Set debug flag for the cluster fitter
799 // Setting fitter parameters
800 fitter->SetInnerWidthFactor(1,1.5);
801 fitter->SetOuterWidthFactor(1,1.5);
802 fitter->SetNmaxOverlaps(5);
804 //fitter->SetChiSqMax(5,kFALSE); //isolated clusters
805 fitter->SetChiSqMax(5,kTRUE); //overlapping clusters
807 Int_t rowrange[2] = {0,AliL3Transform::GetNRows()-1};
809 // Takes input from global hough tracks produced by HT
810 fitter->LoadSeeds(rowrange,kFALSE,iEvent);
814 for(Int_t islice = 0; islice <= 35; islice++)
816 for(Int_t ipatch = 0; ipatch < AliL3Transform::GetNPatches(); ipatch++)
819 fHLTHough->GetMemHandler(ipatch)->Free();
820 fHLTHough->GetMemHandler(ipatch)->Init(islice,ipatch);
821 AliL3DigitRowData *digits = (AliL3DigitRowData *)fHLTHough->GetMemHandler(ipatch)->AliAltroDigits2Memory(ndigits,iEvent);
823 fitter->Init(islice,ipatch);
824 fitter->SetInputData(digits);
825 fitter->FindClusters();
826 fitter->WriteClusters();
830 // Refit of the clusters
832 //The seeds are the input tracks from circle HT
833 AliL3TrackArray *tracks = fitter->GetSeeds();
834 AliL3Fitter *ft = new AliL3Fitter(&vertex,1);
836 ft->LoadClusters("./hlt/fitter/",iEvent,kFALSE);
837 for(Int_t i=0; i<tracks->GetNTracks(); i++)
839 AliL3Track *track = tracks->GetCheckedTrack(i);
841 if(track->GetNHits() < 20) continue;
842 ft->SortTrackClusters(track);
844 track->UpdateToFirstPoint();
848 //Write the final tracks
849 fitter->WriteTracks(20);
853 // remove the event number from the file names
855 sprintf(command, "rename tracks_%d tracks hlt/hough/*.raw", iEvent);
856 gSystem->Exec(command);
857 sprintf(command, "rename tracks_%d tracks hlt/fitter/*.raw", iEvent);
858 gSystem->Exec(command);
859 sprintf(command, "rename points_%d points hlt/fitter/*.raw", iEvent);
860 gSystem->Exec(command);
864 //_____________________________________________________________________________
865 Bool_t AliMonitorProcess::WriteHistos()
867 // write the monitor tree and the monitor histograms to the file
868 // "monitor_<run number>[_<sub_run_number>].root"
869 // if at least fNEventsMin events were monitored
873 // rename tree file and create a new one
880 sprintf(fileName, "monitor_tree_%d.root", fRunNumber);
881 if (fSubRunNumber > 0) {
882 sprintf(fileName, "monitor_tree_%d_%d.root", fRunNumber, fSubRunNumber);
884 if (fNEvents < fNEventsMin) {
885 gSystem->Unlink("monitor_tree.root");
887 gSystem->Rename("monitor_tree.root", fileName);
890 fFile = TFile::Open("monitor_tree.root", "RECREATE");
891 if (!fFile || !fFile->IsOpen()) {
892 AliFatal("could not open file for tree");
894 fTree = new TTree("MonitorTree", "tree for monitoring");
895 for (Int_t iMonitor = 0; iMonitor < fMonitors.GetEntriesFast(); iMonitor++) {
896 ((AliMonitor*) fMonitors[iMonitor])->CreateBranches(fTree);
900 // write the histograms
901 if (fNEvents < fNEventsMin) return kTRUE;
903 if (!fWriteHistoList) {
904 TIterator* iFolder = fTopFolder->GetListOfFolders()->MakeIterator();
905 while (TFolder* folder = (TFolder*) iFolder->Next()) {
906 TIterator* iHisto = folder->GetListOfFolders()->MakeIterator();
907 while (AliMonitorPlot* histo = (AliMonitorPlot*) iHisto->Next()) {
915 Bool_t result = kTRUE;
916 sprintf(fileName, "monitor_%d.root", fRunNumber);
917 if (fSubRunNumber > 0) {
918 sprintf(fileName, "monitor_%d_%d.root", fRunNumber, fSubRunNumber);
920 TFile* file = TFile::Open(fileName, "recreate");
921 if (!file || !file->IsOpen()) {
922 AliError(Form("could not open file %s", fileName));
928 if (file) delete file;
933 //_____________________________________________________________________________
934 void AliMonitorProcess::StartNewRun()
936 // reset the histograms for a new run
938 SetStatus(kResetting);
939 TIterator* iFolder = fTopFolder->GetListOfFolders()->MakeIterator();
940 while (TFolder* folder = (TFolder*) iFolder->Next()) {
941 TIterator* iHisto = folder->GetListOfFolders()->MakeIterator();
942 while (AliMonitorPlot* histo = (AliMonitorPlot*) iHisto->Next()) {
953 //_____________________________________________________________________________
954 void AliMonitorProcess::CheckForConnections()
956 // check if new clients want to connect and add them to the list of sockets
959 while ((socket = fServerSocket->Accept()) != (TSocket*)-1) {
960 socket->SetOption(kNoBlock, 1);
961 char socketType[256];
962 if (socket->Recv(socketType, 255) <= 0) {
963 gSystem->Sleep(1000);
964 if (socket->Recv(socketType, 255) <= 0) {
965 TInetAddress adr = socket->GetInetAddress();
966 AliError(Form("no socket type received - "
967 "disconnect client: %s (%s), port %d",
968 adr.GetHostName(), adr.GetHostAddress(), adr.GetPort()));
973 if (strcmp(socketType, "client") == 0) {
974 fSockets.Add(socket);
975 TInetAddress adr = socket->GetInetAddress();
976 AliInfo(Form("new client: %s (%s), port %d",
977 adr.GetHostName(), adr.GetHostAddress(), adr.GetPort()));
978 if (fNEvents > 0) BroadcastHistos(socket);
979 } else if (strcmp(socketType, "display") == 0) {
980 if (fDisplaySocket) {
981 fDisplaySocket->Close();
982 delete fDisplaySocket;
984 fDisplaySocket = socket;
985 TInetAddress adr = socket->GetInetAddress();
986 AliInfo(Form("new display: %s (%s), port %d",
987 adr.GetHostName(), adr.GetHostAddress(), adr.GetPort()));
989 TInetAddress adr = socket->GetInetAddress();
990 AliError(Form("unknown socket type %s - "
991 "disconnect client: %s (%s), port %d", socketType,
992 adr.GetHostName(), adr.GetHostAddress(), adr.GetPort()));
998 // remove finished or invalid clients
999 for (Int_t iSocket = 0; iSocket < fSockets.GetEntriesFast(); iSocket++) {
1000 socket = (TSocket*) fSockets[iSocket];
1001 if (!socket) continue;
1002 char controlMessage[256];
1003 if (socket->Recv(controlMessage, 255)) {
1004 if (strcmp(controlMessage, "disconnect") == 0) {
1005 TInetAddress adr = socket->GetInetAddress();
1006 AliInfo(Form("disconnect client: %s (%s), port %d",
1007 adr.GetHostName(), adr.GetHostAddress(), adr.GetPort()));
1008 delete fSockets.RemoveAt(iSocket);
1012 if (!socket->IsValid()) {
1013 // remove invalid sockets from the list
1014 TInetAddress adr = socket->GetInetAddress();
1015 AliError(Form("disconnect invalid client: %s (%s), port %d",
1016 adr.GetHostName(), adr.GetHostAddress(), adr.GetPort()));
1017 delete fSockets.RemoveAt(iSocket);
1020 fSockets.Compress();
1023 //_____________________________________________________________________________
1024 void AliMonitorProcess::BroadcastHistos(TSocket* toSocket)
1026 // send the monitor histograms to the clients
1028 SetStatus(kBroadcasting);
1029 TMessage message(kMESS_OBJECT);
1030 message.WriteObject(fTopFolder);
1032 for (Int_t iSocket = 0; iSocket < fSockets.GetEntriesFast(); iSocket++) {
1033 TSocket* socket = (TSocket*) fSockets[iSocket];
1034 if (!socket) continue;
1035 if (toSocket && (socket != toSocket)) continue;
1037 // send control message
1038 if (!socket->IsValid() || (socket->Send("histograms") <= 0)) {
1039 TInetAddress adr = socket->GetInetAddress();
1040 AliError(Form("connection to client failed - "
1041 "disconnect client: %s (%s), port %d",
1042 adr.GetHostName(), adr.GetHostAddress(), adr.GetPort()));
1043 delete fSockets.RemoveAt(iSocket);
1046 // receive control message
1047 char controlMessage[256];
1048 Int_t result = socket->Recv(controlMessage, 255);
1050 gSystem->Sleep(1000); // wait one second and try again
1051 result = socket->Recv(controlMessage, 255);
1054 TInetAddress adr = socket->GetInetAddress();
1055 AliError(Form("no response from client - "
1056 "disconnect client: %s (%s), port %d",
1057 adr.GetHostName(), adr.GetHostAddress(), adr.GetPort()));
1058 delete fSockets.RemoveAt(iSocket);
1061 if (strcmp(controlMessage, "ok") != 0) {
1062 TInetAddress adr = socket->GetInetAddress();
1063 AliError(Form("no \"ok\" message from client - "
1064 "disconnect client: %s (%s), port %d",
1065 adr.GetHostName(), adr.GetHostAddress(), adr.GetPort()));
1066 delete fSockets.RemoveAt(iSocket);
1070 socket->SetOption(kNoBlock, 0);
1071 if (socket->Send(message) < 0) {
1072 // remove the socket from the list if there was an error
1073 TInetAddress adr = socket->GetInetAddress();
1074 AliError(Form("sending histograms failed - "
1075 "disconnect client: %s (%s), port %d",
1076 adr.GetHostName(), adr.GetHostAddress(), adr.GetPort()));
1077 delete fSockets.RemoveAt(iSocket);
1079 gSystem->Sleep(100);
1080 socket->SetOption(kNoBlock, 1);
1083 fSockets.Compress();
1087 //_____________________________________________________________________________
1088 AliMonitorProcess::AliMonitorInterruptHandler::AliMonitorInterruptHandler
1089 (AliMonitorProcess* process):
1090 TSignalHandler(kSigUser1, kFALSE),
1093 // constructor: set process
1096 //_____________________________________________________________________________
1097 AliMonitorProcess::AliMonitorInterruptHandler::AliMonitorInterruptHandler
1098 (const AliMonitorInterruptHandler& handler):
1099 TSignalHandler(handler)
1103 AliFatal("copy constructor not implemented");
1106 //_____________________________________________________________________________
1107 AliMonitorProcess::AliMonitorInterruptHandler&
1108 AliMonitorProcess::AliMonitorInterruptHandler::operator =
1109 (const AliMonitorInterruptHandler& /*handler*/)
1111 // assignment operator
1113 AliFatal("assignment operator not implemented");
1117 //_____________________________________________________________________________
1118 Bool_t AliMonitorProcess::AliMonitorInterruptHandler::Notify()
1120 // interrupt signal -> stop process
1122 AliInfo("the monitoring process will be stopped.");