]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/MUON/macros/MonitorRawData.C
Making changes so this macro works better in the ACR.
[u/mrichter/AliRoot.git] / HLT / MUON / macros / MonitorRawData.C
1 /**************************************************************************
2  * This file is property of and copyright by the ALICE HLT Project        *
3  * All rights reserved.                                                   *
4  *                                                                        *
5  * Primary Authors:                                                       *
6  *   Artur Szostak <artursz@iafrica.com>                                  *
7  *                                                                        *
8  * Permission to use, copy, modify and distribute this software and its   *
9  * documentation strictly for non-commercial purposes is hereby granted   *
10  * without fee, provided that the above copyright notice appears in all   *
11  * copies and that both the copyright notice and this permission notice   *
12  * appear in the supporting documentation. The authors make no claims     *
13  * about the suitability of this software for any purpose. It is          *
14  * provided "as is" without express or implied warranty.                  *
15  **************************************************************************/
16
17 #if !defined(__CINT__) || defined(__MAKECINT__)
18 #define USE_ROOT
19 #include "AliHLTHOMERReader.h"
20 #include "AliHLTMessage.h"
21 #include "Riostream.h"
22 #include "Rtypes.h"
23 #include "TH1D.h"
24 #include "TObject.h"
25 #include "TClass.h"
26 #include "TDirectory.h"
27 #include "TTimeStamp.h"
28 #include "TSystem.h"
29 #include "TCanvas.h"
30 #include <errno.h>
31 #else
32 #error You must compile this macro. Try the following command "> aliroot MonitorRawData.C++"
33 #endif
34
35
36 TH1D* gErrorHist = NULL;
37 TH1D* gManuHist = NULL;
38 TH1D* gSignalHist = NULL;
39
40
41 const char* DataTypeToString(homer_uint64 type)
42 {
43         union
44         {
45                 homer_uint64 val;
46                 char bytes[8];
47         };
48         val = type;
49         
50         static char str[9];
51         for (int i = 0; i < 8; i++)
52         {
53                 str[i] = bytes[7-i];
54         }
55         str[8] = '\0'; // Null terminate the string.
56         return str;
57 }
58
59
60 const char* OriginToString(homer_uint32 origin)
61 {
62         union
63         {
64                 homer_uint32 val;
65                 char bytes[4];
66         };
67         val = origin;
68         
69         static char str[5];
70         for (int i = 0; i < 4; i++)
71         {
72                 str[i] = bytes[3-i];
73         }
74         str[4] = '\0'; // Null terminate the string.
75         return str;
76 }
77
78
79 bool UpdateHists(AliHLTHOMERReader& homerReader, const char* hostname, UShort_t port, bool addHists)
80 {
81         /// This routine just reads all the data blocks from the HOMER reader interface,
82         /// and checks if they are histogram objects we can deal with.
83         /// Any histogram objects with the following strings in their names:
84         /// "rawDataErrors", "manuDistrib" or "signalDistrib"
85         /// will be added to the global cumulative histograms.
86
87         int result = homerReader.ReadNextEvent(3000000);  // 3 second timeout.
88         TTimeStamp now;
89         if (result != 0)
90         {
91                 if (result != ETIMEDOUT)
92                 {
93                         cerr << "ERROR: Could not read another event from HLT on "
94                                 << hostname << ":" << port << "\t" << now.AsString() << endl;
95                         return false;
96                 }
97                 else
98                 {
99                         cerr << "Timed out when trying to read from HLT on "
100                                 << hostname << ":" << port << "\t" << now.AsString() << endl;
101                         return true;
102                 }
103         }       
104         
105         cout << "Reading stats up to event: " << homerReader.GetEventID()
106                 << "\t" << now.AsString() << endl;
107         
108         for (unsigned long n = 0; n < homerReader.GetBlockCnt(); n++)
109         {
110                 char* buffer = new char[homerReader.GetBlockDataLength(n)];
111                 memcpy(buffer, homerReader.GetBlockData(n), homerReader.GetBlockDataLength(n));
112                 AliHLTMessage msg(buffer, homerReader.GetBlockDataLength(n));
113                 TClass* objclass = msg.GetClass();
114                 if (objclass == NULL)
115                 {
116                         cerr << "WARNING: Do not know how to handle block " << n
117                                 << ", block type = " << DataTypeToString(homerReader.GetBlockDataType(n))
118                                 << ", origin = " << OriginToString(homerReader.GetBlockDataOrigin(n))
119                                 << ", but unknown class type. Skipping block." << endl;
120                         //delete buffer;
121                         continue;
122                 }
123                 TObject* obj = msg.ReadObject(objclass);
124                 if (obj == NULL)
125                 {
126                         cerr << "WARNING: Could not read object from block " << n
127                                 << ", class type = " << objclass->GetName()
128                                 << ". Skipping block." << endl;
129                         delete buffer;
130                         continue;
131                 }
132                 if (obj->IsA() != TH1D::Class())
133                 {
134                         cerr << "WARNING: Do not know how to handle object of type "
135                                 << obj->ClassName() << " recevied in block " << n
136                                 << ". Skipping block." << endl;
137                         delete buffer;
138                         continue;
139                 }
140                 
141                 TH1D* hist = static_cast<TH1D*>(obj);
142                 TString name = hist->GetName();
143                 if (name.Contains("rawDataErrors"))
144                 {
145                         if (! addHists) gErrorHist->Reset("M");
146                         gErrorHist->Add(hist);
147                 }
148                 else if (name.Contains("manuDistrib"))
149                 {
150                         if (! addHists) gManuHist->Reset("M");
151                         gManuHist->Add(hist);
152                 }
153                 else if (name.Contains("signalDistrib"))
154                 {
155                         if (! addHists) gSignalHist->Reset("M");
156                         gSignalHist->Add(hist);
157                 }
158                 else
159                 {
160                         cerr << "WARNING: Do not know how to handle histogram " << name
161                                 << " found in data block " << n << endl;
162                 }
163                 delete buffer;
164         }
165         
166         return true;
167 }
168
169
170 void MonitorRawData(const char* hostname = "alihlt-dcs0", UShort_t port = 58784, bool addHists = false)
171 {
172         if (gErrorHist == NULL && gDirectory->Get("errorHist") == NULL)
173         {
174                 gErrorHist = new TH1D("errorHist", "Error codes found in raw data", 40, 0.5, 40.5);
175                 gErrorHist->SetXTitle("Error code");
176                 gErrorHist->SetYTitle("Number of errors");
177         }
178         if (gManuHist == NULL && gDirectory->Get("manuHist") == NULL)
179         {
180                 gManuHist = new TH1D("manuHist", "Distribution of signals found per MANU", 2048, -0.5, 2047.5);
181                 gManuHist->SetXTitle("MANU number (as seen in raw data)");
182                 gManuHist->SetYTitle("Number of signals received.");
183         }
184         if (gSignalHist == NULL && gDirectory->Get("signalHist") == NULL)
185         {
186                 gSignalHist = new TH1D("signalHist", "Distribution of ADC signal values", 4096, -0.5, 4095.5);
187                 gSignalHist->SetXTitle("Channels");
188                 gSignalHist->SetYTitle("dN/dChannel");
189         }
190         
191         TCanvas* c1 = new TCanvas("errorCanvas", "Error histogram", 0, 0, 600, 450);
192         gErrorHist->Draw();
193         TCanvas* c2 = new TCanvas("manuCanvas", "MANU histogram", 610, 00, 600, 450);
194         if (gManuHist->GetEntries() > 0) c2->SetLogy();
195         gManuHist->Draw();
196         TCanvas* c3 = new TCanvas("signalCanvas", "Signal histogram", 0, 480, 600, 450);
197         if (gSignalHist->GetEntries() > 0) c3->SetLogy();
198         gSignalHist->Draw();
199         
200         for (;;)
201         {
202                 bool updateOk = true;
203
204                 AliHLTHOMERReader homerReader(hostname, port);
205                 int status = homerReader.GetConnectionStatus();
206                 if (status == 0)
207                 {
208                         try
209                         {
210                                 updateOk = UpdateHists(homerReader, hostname, port, addHists);
211                         }
212                         catch (...)
213                         {
214                                 cerr << "ERROR: exception occured. Trying to recover and continue..." << endl;
215                         }
216                 }
217                 else
218                 {
219                         cerr << "ERROR: Could not connect to HLT on " << hostname << ":" << port << endl;
220                 }
221
222                 // Clear histograms on errors. This is usualy due to end of run.
223                 if (! updateOk || status != 0)
224                 {
225                         gErrorHist->Reset("M");
226                         gManuHist->Reset("M");
227                         gSignalHist->Reset("M");
228                         c1->cd();
229                         gErrorHist->Draw();
230                         c2->cd();
231                         c2->SetLogy(kFALSE);
232                         gManuHist->Draw();
233                         c3->cd();
234                         c3->SetLogy(kFALSE);
235                         gSignalHist->Draw();
236                 }
237
238                 if (gManuHist->GetEntries() > 0) c2->SetLogy();
239                 if (gSignalHist->GetEntries() > 0) c3->SetLogy();
240
241                 c1->Update();
242                 c2->Update();
243                 c3->Update();
244
245                 // Effectively wait for 2 seconds but dispatch ROOT events while waiting.
246                 for (int i = 0; i < 200; i++)
247                 {
248                         gSystem->DispatchOneEvent(kTRUE);
249                         gSystem->Sleep(10);
250                 }
251         }
252 }
253