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 "AliRawReaderRoot.h"
35 #include "AliLoader.h"
38 #include "AliTPCclustererMI.h"
39 #include "AliTPCtrackerMI.h"
41 #include "AliITSclustererV2.h"
42 #include "AliITStrackerV2.h"
43 #include "AliITSLoader.h"
44 #include "AliV0vertexer.h"
47 #include <TServerSocket.h>
49 #include <TGridResult.h>
52 #include <AliLevel3.h>
53 #include <AliL3Transform.h>
56 ClassImp(AliMonitorProcess)
59 const Int_t AliMonitorProcess::fgkPort = 9327;
62 //_____________________________________________________________________________
63 AliMonitorProcess::AliMonitorProcess(const char* alienDir,
64 const char* fileNameGalice)
66 // initialize the monitoring process and the monitor histograms
68 fGrid = TGrid::Connect("alien", gSystem->Getenv("USER"));
69 if (!fGrid || fGrid->IsZombie() || !fGrid->IsConnected()) {
71 Fatal("AliMonitorProcess", "could not connect to alien");
74 fLogicalFileName = "";
77 fRunLoader = AliRunLoader::Open(fileNameGalice);
78 if (!fRunLoader) Fatal("AliMonitorProcess",
79 "could not get run loader from file %s",
82 fRunLoader->CdGAFile();
83 fTPCParam = AliTPC::LoadTPCParam(gFile);
84 if (!fTPCParam) Fatal("AliMonitorProcess", "could not load TPC parameters");
86 fRunLoader->LoadgAlice();
87 gAlice = fRunLoader->GetAliRun();
88 if (!gAlice) Fatal("AliMonitorProcess", "no gAlice object found");
89 AliITS* its = (AliITS*) gAlice->GetModule("ITS");
90 if (!its) Fatal("AliMonitorProcess", "no ITS detector found");
91 fITSgeom = its->GetITSgeom();
92 if (!fITSgeom) Fatal("AliMonitorProcess", "could not load ITS geometry");
95 // Init TPC parameters for HLT
96 Bool_t isinit=AliL3Transform::Init(const_cast<char*>(fileNameGalice),kTRUE);
98 cerr << "Could not create transform settings, please check log for error messages!" << endl;
107 fWriteHistoList = kFALSE;
109 fTopFolder = new TFolder("Monitor", "monitor histograms");
110 fTopFolder->SetOwner(kTRUE);
112 fMonitors.Add(new AliMonitorTPC(fTPCParam));
113 fMonitors.Add(new AliMonitorITS(fITSgeom));
114 fMonitors.Add(new AliMonitorV0s);
116 fMonitors.Add(new AliMonitorHLT(fTPCParam));
119 for (Int_t iMonitor = 0; iMonitor < fMonitors.GetEntriesFast(); iMonitor++) {
120 ((AliMonitor*) fMonitors[iMonitor])->CreateHistos(fTopFolder);
123 TIterator* iFolder = fTopFolder->GetListOfFolders()->MakeIterator();
124 while (TFolder* folder = (TFolder*) iFolder->Next()) folder->SetOwner(kTRUE);
127 fFile = TFile::Open("monitor_tree.root", "RECREATE");
128 if (!fFile || !fFile->IsOpen()) {
129 Fatal("AliMonitorProcess", "could not open file for tree");
131 fTree = new TTree("MonitorTree", "tree for monitoring");
132 for (Int_t iMonitor = 0; iMonitor < fMonitors.GetEntriesFast(); iMonitor++) {
133 ((AliMonitor*) fMonitors[iMonitor])->CreateBranches(fTree);
137 fServerSocket = new TServerSocket(fgkPort, kTRUE);
138 fServerSocket->SetOption(kNoBlock, 1);
139 fDisplaySocket = NULL;
140 CheckForConnections();
149 //_____________________________________________________________________________
150 AliMonitorProcess::AliMonitorProcess(const AliMonitorProcess& process) :
153 Fatal("AliMonitorProcess", "copy constructor not implemented");
156 //_____________________________________________________________________________
157 AliMonitorProcess& AliMonitorProcess::operator = (const AliMonitorProcess&
160 Fatal("operator =", "assignment operator not implemented");
164 //_____________________________________________________________________________
165 AliMonitorProcess::~AliMonitorProcess()
173 delete fServerSocket;
175 delete fDisplaySocket;
182 gSystem->Unlink("monitor_tree.root");
190 //_____________________________________________________________________________
191 const char* AliMonitorProcess::GetRevision()
197 //_____________________________________________________________________________
198 void AliMonitorProcess::SetStatus(EStatus status)
200 // set the current status and process system events
203 gSystem->ProcessEvents();
207 //_____________________________________________________________________________
208 void AliMonitorProcess::Run()
210 // run the monitor process:
211 // check for a raw data file, process the raw data file and delete it
217 while (!CheckForNewFile()) {
218 CheckForConnections();
220 if (fStopping) break;
223 if (fStopping) break;
235 //_____________________________________________________________________________
236 void AliMonitorProcess::Stop()
238 // set the fStopping flag to terminate the monitor process after the current
239 // event was processed
241 if (GetStatus() != kStopped) fStopping = kTRUE;
245 //_____________________________________________________________________________
246 void AliMonitorProcess::ProcessFile(const char* fileName)
248 // create a file with monitor histograms for a single file
250 if (GetStatus() != kStopped) {
251 Error("ProcessFile", "ProcessFile can not be called"
252 " while the monitor process is running");
256 fFileName = fileName;
257 Int_t nEventMin = fNEventsMin;
261 fNEventsMin = nEventMin;
266 //_____________________________________________________________________________
267 Bool_t AliMonitorProcess::CheckForNewFile()
269 // check whether a new file was registered in alien
271 TGridResult* result = fGrid->Ls();
276 while (const char* entry = result->Next()) {
277 // entry = host_date_time.root
278 TString entryCopy(entry);
279 char* p = const_cast<char*>(entryCopy.Data());
280 if (!strtok(p, "_") || !p) continue; // host name
281 char* dateStr = strtok(NULL, "_");
282 if (!dateStr || !p) continue;
283 char* timeStr = strtok(NULL, ".");
284 if (!timeStr || !p) continue;
285 Long_t date = atoi(dateStr);
286 Long_t time = atoi(timeStr);
288 if ((date > maxDate) || ((date == maxDate) && (time > maxTime))) {
295 if (maxDate < 0) return kFALSE; // no files found
296 if (fLogicalFileName.CompareTo(fileName) == 0) return kFALSE; // no new file
298 fLogicalFileName = fileName;
299 TGridResult* result2 = fGrid->GetPhysicalFileNames(fLogicalFileName.Data());
300 fFileName = result2->Next();
305 //_____________________________________________________________________________
306 Bool_t AliMonitorProcess::ProcessFile()
308 // loop over all events in the raw data file, run the reconstruction
309 // and fill the monitor histograms
311 Int_t nEvents = GetNumberOfEvents(fFileName);
312 if (nEvents <= 0) return kFALSE;
313 Info("ProcessFile", "found %d event(s) in file %s",
314 nEvents, fFileName.Data());
316 CreateHLT(fFileName);
319 // loop over the events
320 for (Int_t iEvent = 0; iEvent < nEvents; iEvent++) {
322 fRunLoader->SetEventNumber(0);
323 AliRawReaderRoot rawReader(fFileName, iEvent);
324 if (fStopping) break;
325 if (rawReader.GetRunNumber() != fRunNumber) {
328 fRunNumber = rawReader.GetRunNumber();
329 fEventNumber[0] = rawReader.GetEventId()[0];
330 fEventNumber[1] = rawReader.GetEventId()[1];
332 if (fStopping) break;
335 if (!ReconstructTPC(&rawReader)) return kFALSE;
336 if (fStopping) break;
337 if (!ReconstructITS(&rawReader)) return kFALSE;
338 if (fStopping) break;
339 if (!ReconstructV0s()) return kFALSE;
340 if (fStopping) break;
341 if (!ReconstructHLT(iEvent)) return kFALSE;
342 if (fStopping) break;
344 if (fDisplaySocket) fDisplaySocket->Send("new event");
346 Info("ProcessFile", "filling histograms...");
348 for (Int_t iMonitor = 0; iMonitor < fMonitors.GetEntriesFast(); iMonitor++) {
349 ((AliMonitor*) fMonitors[iMonitor])->FillHistos(fRunLoader, &rawReader);
350 if (fStopping) break;
352 if (fStopping) break;
354 Info("ProcessFile", "updating histograms...");
355 SetStatus(kUpdating);
356 TIterator* iFolder = fTopFolder->GetListOfFolders()->MakeIterator();
357 while (TFolder* folder = (TFolder*) iFolder->Next()) {
358 TIterator* iHisto = folder->GetListOfFolders()->MakeIterator();
359 while (AliMonitorPlot* histo = (AliMonitorPlot*) iHisto->Next()) {
365 if (fStopping) break;
367 Info("ProcessFile", "filling the tree...");
370 Info("ProcessFile", "broadcasting histograms...");
371 CheckForConnections();
375 if (fStopping) break;
386 //_____________________________________________________________________________
387 void AliMonitorProcess::Reset()
389 // write the current histograms to a file and reset them
391 if (fSubRunNumber == 0) fSubRunNumber++;
398 //_____________________________________________________________________________
399 UInt_t AliMonitorProcess::GetEventPeriodNumber() const
401 // get the period number from the event id
403 return (fEventNumber[1] >> 4);
406 //_____________________________________________________________________________
407 UInt_t AliMonitorProcess::GetEventOrbitNumber() const
409 // get the orbit number from the event id
411 return ((fEventNumber[1] & 0x000F) << 20) + (fEventNumber[0] >> 12);
414 //_____________________________________________________________________________
415 UInt_t AliMonitorProcess::GetEventBunchNumber() const
417 // get the bunch number from the event id
419 return (fEventNumber[0] % 0x0FFF);
422 //_____________________________________________________________________________
423 Int_t AliMonitorProcess::GetNumberOfEvents(const char* fileName) const
425 // determine the number of events in the given raw data file
429 TFile* file = TFile::Open(fileName);
430 if (!file || !file->IsOpen()) {
431 Error("GetNumberOfEvents", "could not open file %s", fileName);
432 if (file) delete file;
436 TTree* tree = (TTree*) file->Get("RAW");
438 Error("GetNumberOfEvents", "could not find tree with raw data");
440 nEvents = (Int_t) tree->GetEntries();
448 //_____________________________________________________________________________
449 Bool_t AliMonitorProcess::ReconstructTPC(AliRawReader* rawReader)
451 // find TPC clusters and tracks
455 AliLoader* tpcLoader = fRunLoader->GetLoader("TPCLoader");
457 Error("ReconstructTPC", "no TPC loader found");
460 gSystem->Unlink("TPC.RecPoints.root");
461 gSystem->Unlink("TPC.Tracks.root");
464 Info("ReconstructTPC", "reconstructing clusters...");
465 tpcLoader->LoadRecPoints("recreate");
466 AliTPCclustererMI clusterer(fTPCParam);
467 tpcLoader->MakeRecPointsContainer();
468 clusterer.SetOutput(tpcLoader->TreeR());
469 clusterer.Digits2Clusters(rawReader);
470 tpcLoader->WriteRecPoints("OVERWRITE");
473 Info("ReconstructTPC", "reconstructing tracks...");
474 tpcLoader->LoadTracks("recreate");
476 AliTPCtrackerMI tracker(fTPCParam);
477 tracker.Clusters2Tracks();
480 tpcLoader->UnloadRecPoints();
481 tpcLoader->UnloadTracks();
485 //_____________________________________________________________________________
486 Bool_t AliMonitorProcess::ReconstructITS(AliRawReader* rawReader)
488 // find ITS clusters and tracks
492 AliLoader* itsLoader = fRunLoader->GetLoader("ITSLoader");
494 Error("ReconstructITS", "no ITS loader found");
497 AliLoader* tpcLoader = fRunLoader->GetLoader("TPCLoader");
499 Error("ReconstructITS", "no TPC loader found");
502 gSystem->Unlink("ITS.RecPoints.root");
503 gSystem->Unlink("ITS.Tracks.root");
506 Info("ReconstructITS", "reconstructing clusters...");
507 itsLoader->LoadRecPoints("recreate");
508 AliITSclustererV2 clusterer(fITSgeom);
509 itsLoader->MakeRecPointsContainer();
510 clusterer.Digits2Clusters(rawReader);
513 Info("ReconstructITS", "reconstructing tracks...");
514 itsLoader->LoadTracks("recreate");
515 itsLoader->MakeTracksContainer();
516 tpcLoader->LoadTracks();
517 AliITStrackerV2 tracker(fITSgeom);
518 tracker.LoadClusters(itsLoader->TreeR());
519 tracker.Clusters2Tracks(tpcLoader->TreeT(), itsLoader->TreeT());
520 tracker.UnloadClusters();
521 itsLoader->WriteTracks("OVERWRITE");
523 itsLoader->UnloadRecPoints();
524 itsLoader->UnloadTracks();
525 tpcLoader->UnloadTracks();
529 //_____________________________________________________________________________
530 Bool_t AliMonitorProcess::ReconstructV0s()
536 AliITSLoader* itsLoader = (AliITSLoader*) fRunLoader->GetLoader("ITSLoader");
538 Error("ReconstructV0", "no ITS loader found");
541 gSystem->Unlink("ITS.V0s.root");
544 Info("ReconstructV0s", "reconstructing V0s...");
545 itsLoader->LoadTracks("read");
546 itsLoader->LoadV0s("recreate");
547 AliV0vertexer vertexer;
548 TTree* tracks = itsLoader->TreeT();
550 Error("ReconstructV0s", "no ITS tracks tree found");
553 if (!itsLoader->TreeV0()) itsLoader->MakeTree("V0");
554 TTree* v0s = itsLoader->TreeV0();
555 vertexer.Tracks2V0vertices(tracks, v0s);
556 itsLoader->WriteV0s("OVERWRITE");
558 itsLoader->UnloadTracks();
559 itsLoader->UnloadV0s();
563 //_____________________________________________________________________________
565 void AliMonitorProcess::CreateHLT(const char* fileName)
568 // create the HLT (Level3) object
570 if (fHLT) delete fHLT;
573 strcpy(name, fileName);
574 fHLT = new AliLevel3(name);
575 fHLT->Init("./", AliLevel3::kRaw, 1);
577 fHLT->SetClusterFinderParam(0, 0, kTRUE);
579 Int_t phiSegments = 50;
580 Int_t etaSegments = 100;
581 Int_t trackletlength = 3;
582 Int_t tracklength = 40;//40 or 5
583 Int_t rowscopetracklet = 2;
584 Int_t rowscopetrack = 2;
585 Double_t minPtFit = 0;
586 Double_t maxangle = 1.31;
587 Double_t goodDist = 5;
588 Double_t maxphi = 100;
589 Double_t maxeta = 100;
590 Double_t hitChi2Cut = 15;//100 or 15
591 Double_t goodHitChi2 = 5;//20 or 5
592 Double_t trackChi2Cut = 10;//50 or 10
593 fHLT->SetTrackerParam(phiSegments, etaSegments,
594 trackletlength, tracklength,
595 rowscopetracklet, rowscopetrack,
596 minPtFit, maxangle, goodDist, hitChi2Cut,
597 goodHitChi2, trackChi2Cut, 50, maxphi, maxeta, kTRUE);
599 fHLT->WriteFiles("./hlt/");
603 //_____________________________________________________________________________
604 Bool_t AliMonitorProcess::ReconstructHLT(
612 // run the HLT cluster and track finder
617 Warning("ReconstructHLT", "the code was compiled without HLT support");
621 gSystem->Exec("rm -rf hlt");
622 gSystem->MakeDirectory("hlt");
623 if (!fHLT) return kFALSE;
625 fHLT->ProcessEvent(0, 35, iEvent);
627 // remove the event number from the file names
629 sprintf(command, "rename points_%d points hlt/*.raw", iEvent);
630 gSystem->Exec(command);
631 sprintf(command, "rename tracks_tr_%d tracks_tr hlt/*.raw", iEvent);
632 gSystem->Exec(command);
633 sprintf(command, "rename tracks_gl_%d tracks_gl hlt/*.raw", iEvent);
634 gSystem->Exec(command);
635 sprintf(command, "rename tracks_%d tracks hlt/*.raw", iEvent);
636 gSystem->Exec(command);
642 //_____________________________________________________________________________
643 Bool_t AliMonitorProcess::WriteHistos()
645 // write the monitor tree and the monitor histograms to the file
646 // "monitor_<run number>[_<sub_run_number>].root"
647 // if at least fNEventsMin events were monitored
651 // rename tree file and create a new one
658 sprintf(fileName, "monitor_tree_%d.root", fRunNumber);
659 if (fSubRunNumber > 0) {
660 sprintf(fileName, "monitor_tree_%d_%d.root", fRunNumber, fSubRunNumber);
662 if (fNEvents < fNEventsMin) {
663 gSystem->Unlink("monitor_tree.root");
665 gSystem->Rename("monitor_tree.root", fileName);
668 fFile = TFile::Open("monitor_tree.root", "RECREATE");
669 if (!fFile || !fFile->IsOpen()) {
670 Fatal("WriteHistos", "could not open file for tree");
672 fTree = new TTree("MonitorTree", "tree for monitoring");
673 for (Int_t iMonitor = 0; iMonitor < fMonitors.GetEntriesFast(); iMonitor++) {
674 ((AliMonitor*) fMonitors[iMonitor])->CreateBranches(fTree);
678 // write the histograms
679 if (fNEvents < fNEventsMin) return kTRUE;
681 if (!fWriteHistoList) {
682 TIterator* iFolder = fTopFolder->GetListOfFolders()->MakeIterator();
683 while (TFolder* folder = (TFolder*) iFolder->Next()) {
684 TIterator* iHisto = folder->GetListOfFolders()->MakeIterator();
685 while (AliMonitorPlot* histo = (AliMonitorPlot*) iHisto->Next()) {
693 Bool_t result = kTRUE;
694 sprintf(fileName, "monitor_%d.root", fRunNumber);
695 if (fSubRunNumber > 0) {
696 sprintf(fileName, "monitor_%d_%d.root", fRunNumber, fSubRunNumber);
698 TFile* file = TFile::Open(fileName, "recreate");
699 if (!file || !file->IsOpen()) {
700 Error("WriteHistos", "could not open file %s", fileName);
706 if (file) delete file;
711 //_____________________________________________________________________________
712 void AliMonitorProcess::StartNewRun()
714 // reset the histograms for a new run
716 SetStatus(kResetting);
717 TIterator* iFolder = fTopFolder->GetListOfFolders()->MakeIterator();
718 while (TFolder* folder = (TFolder*) iFolder->Next()) {
719 TIterator* iHisto = folder->GetListOfFolders()->MakeIterator();
720 while (AliMonitorPlot* histo = (AliMonitorPlot*) iHisto->Next()) {
731 //_____________________________________________________________________________
732 void AliMonitorProcess::CheckForConnections()
734 // check if new clients want to connect and add them to the list of sockets
736 TMessage message(kMESS_OBJECT);
737 message.WriteObject(fTopFolder);
738 SetStatus(kConnecting);
741 while ((socket = fServerSocket->Accept()) != (TSocket*)-1) {
742 socket->SetOption(kNoBlock, 1);
743 char socketType[256];
744 if (!socket->Recv(socketType, 255)) continue;
745 if (strcmp(socketType, "client") == 0) {
746 if ((fNEvents == 0) || (socket->Send(message) >= 0)) {
747 fSockets.Add(socket);
748 TInetAddress adr = socket->GetInetAddress();
749 Info("CheckForConnections", "new client:\n %s (%s), port %d\n",
750 adr.GetHostName(), adr.GetHostAddress(), adr.GetPort());
752 } else if (strcmp(socketType, "display") == 0) {
753 if (fDisplaySocket) {
754 fDisplaySocket->Close();
755 delete fDisplaySocket;
757 fDisplaySocket = socket;
758 fDisplaySocket->SetOption(kNoBlock, 1);
759 TInetAddress adr = socket->GetInetAddress();
760 Info("CheckForConnections", "new display:\n %s (%s), port %d\n",
761 adr.GetHostName(), adr.GetHostAddress(), adr.GetPort());
765 for (Int_t iSocket = 0; iSocket < fSockets.GetEntriesFast(); iSocket++) {
766 socket = (TSocket*) fSockets[iSocket];
767 if (!socket) continue;
768 // remove finished client
770 if (socket->Recv(str, 255)) {
771 TString socketMessage(str);
772 if(socketMessage.CompareTo("Finished") == 0) {
773 TInetAddress adr = socket->GetInetAddress();
774 Info("CheckForConnections",
775 "disconnect finished client:\n %s (%s), port %d\n",
776 adr.GetHostName(), adr.GetHostAddress(), adr.GetPort());
777 delete fSockets.RemoveAt(iSocket);
781 if (!socket->IsValid()) {
782 // remove invalid sockets from the list
783 TInetAddress adr = socket->GetInetAddress();
784 Info("BroadcastHistos", "disconnect client:\n %s (%s), port %d\n",
785 adr.GetHostName(), adr.GetHostAddress(), adr.GetPort());
786 delete fSockets.RemoveAt(iSocket);
792 //_____________________________________________________________________________
793 void AliMonitorProcess::BroadcastHistos()
795 // send the monitor histograms to the clients
797 SetStatus(kBroadcasting);
798 TMessage message(kMESS_OBJECT);
799 message.WriteObject(fTopFolder);
801 for (Int_t iSocket = 0; iSocket < fSockets.GetEntriesFast(); iSocket++) {
802 TSocket* socket = (TSocket*) fSockets[iSocket];
803 if (!socket) continue;
804 socket->SetOption(kNoBlock, 0);
805 if (!socket->IsValid() || (socket->Send(message) < 0)) {
806 // remove the socket from the list if there was an error
807 TInetAddress adr = socket->GetInetAddress();
808 Info("BroadcastHistos", "disconnect client:\n %s (%s), port %d\n",
809 adr.GetHostName(), adr.GetHostAddress(), adr.GetPort());
810 delete fSockets.RemoveAt(iSocket);
812 socket->SetOption(kNoBlock, 1);