]>
Commit | Line | Data |
---|---|---|
a63da6d6 | 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 | ||
1d8ae082 | 17 | // $Id: $ |
a63da6d6 | 18 | |
19 | /// | |
20 | /// @file AliHLTMUONRawDataHistoComponent.cxx | |
21 | /// @author Artur Szostak <artursz@iafrica.com> | |
22 | /// @date 30 April 2008 | |
23 | /// @brief Implementation of the raw data histogramming component for dHLT. | |
24 | /// | |
1d8ae082 | 25 | /// The class implements a component for checking basic statistics and errors in |
26 | /// raw data from the muon spectrometer. It is useful for basic monitoring. | |
a63da6d6 | 27 | |
28 | #include "AliHLTMUONRawDataHistoComponent.h" | |
29 | #include "AliHLTMUONConstants.h" | |
30 | #include "AliHLTMUONUtils.h" | |
31 | #include "AliHLTDataTypes.h" | |
32 | #include "AliCDBEntry.h" | |
33 | #include "AliCDBManager.h" | |
34 | #include "AliRawDataHeader.h" | |
a63da6d6 | 35 | #include "TTimeStamp.h" |
36 | #include <cstdlib> | |
37 | #include <cstring> | |
38 | #include <cerrno> | |
39 | #include <cmath> | |
40 | #include <new> | |
41 | ||
42 | ||
43 | // Helper type for memory allocation. | |
44 | typedef const AliHLTMUONMansoTrackStruct* AliHLTMUONMansoTrackStructP; | |
45 | ||
46 | ||
47 | ClassImp(AliHLTMUONRawDataHistoComponent); | |
48 | ||
49 | ||
50 | AliHLTMUONRawDataHistoComponent::AliHLTMUONRawDataHistoComponent() : | |
51 | AliHLTMUONProcessor(), | |
7989fd8e | 52 | fTrackerDecoder(), |
53 | fTriggerDecoder(), | |
a63da6d6 | 54 | fLastPublishTime(-1), |
55 | fCurrentEventTime(-1), | |
a35ef960 | 56 | fPublishDelay(0), |
7989fd8e | 57 | fSuppressEmptyHists(false), |
a35ef960 | 58 | fProcessDataEventsOnly(false), |
59 | fClearAfterPublish(false) | |
a63da6d6 | 60 | { |
61 | /// Default constructor initialises all histogram object pointers to NULL. | |
62 | ||
63 | for (int i = 0; i < 22; i++) | |
64 | { | |
65 | fErrorHist[i] = NULL; | |
66 | } | |
67 | for (int i = 0; i < 20; i++) | |
68 | { | |
69 | fManuHist[i] = NULL; | |
70 | fSignalHist[i] = NULL; | |
71 | } | |
7989fd8e | 72 | |
73 | fTrackerDecoder.ExitOnError(false); | |
74 | fTrackerDecoder.TryRecover(false); | |
75 | fTrackerDecoder.SendDataOnParityError(true); | |
76 | fTrackerDecoder.AutoDetectTrailer(true); | |
77 | fTrackerDecoder.CheckForTrailer(true); | |
78 | ||
79 | fTriggerDecoder.ExitOnError(false); | |
80 | fTriggerDecoder.TryRecover(false); | |
81 | fTriggerDecoder.AutoDetectScalars(false); | |
a63da6d6 | 82 | } |
83 | ||
84 | ||
85 | AliHLTMUONRawDataHistoComponent::~AliHLTMUONRawDataHistoComponent() | |
86 | { | |
87 | /// Default destructor deletes all histogram objects if they are still allocated. | |
88 | ||
89 | } | |
90 | ||
91 | ||
92 | const char* AliHLTMUONRawDataHistoComponent::GetComponentID() | |
93 | { | |
a63da6d6 | 94 | /// Inherited from AliHLTComponent. Returns the component ID. |
a63da6d6 | 95 | |
96 | return AliHLTMUONConstants::RawDataHistogrammerId(); | |
97 | } | |
98 | ||
99 | ||
100 | void AliHLTMUONRawDataHistoComponent::GetInputDataTypes(AliHLTComponentDataTypeList& list) | |
101 | { | |
a63da6d6 | 102 | /// Inherited from AliHLTProcessor. Returns the list of expected input data types. |
a63da6d6 | 103 | |
104 | assert( list.empty() ); | |
105 | list.push_back( AliHLTMUONConstants::DDLRawDataType() ); | |
106 | } | |
107 | ||
108 | ||
109 | AliHLTComponentDataType AliHLTMUONRawDataHistoComponent::GetOutputDataType() | |
110 | { | |
111 | /// Inherited from AliHLTComponent. Returns kAliHLTHistogramDataTypeID. | |
112 | ||
113 | return AliHLTMUONConstants::HistogramDataType(); | |
114 | } | |
115 | ||
116 | ||
117 | void AliHLTMUONRawDataHistoComponent::GetOutputDataSize( | |
118 | unsigned long& constBase, double& inputMultiplier | |
119 | ) | |
120 | { | |
a63da6d6 | 121 | /// Inherited from AliHLTComponent. Returns an estimate of the expected output data size. |
a63da6d6 | 122 | |
95c6f20f | 123 | constBase = (sizeof(TH1D)+50*sizeof(double))*22 + (sizeof(TH1D)+1024*4*sizeof(double))*20*2; |
a63da6d6 | 124 | inputMultiplier = 0; |
125 | } | |
126 | ||
127 | ||
128 | AliHLTComponent* AliHLTMUONRawDataHistoComponent::Spawn() | |
129 | { | |
a63da6d6 | 130 | /// Inherited from AliHLTComponent. Creates a new object instance. |
a63da6d6 | 131 | |
132 | return new AliHLTMUONRawDataHistoComponent; | |
133 | } | |
134 | ||
135 | ||
136 | bool AliHLTMUONRawDataHistoComponent::IgnoreArgument(const char* arg) const | |
137 | { | |
138 | /// Return true if the argument is one of -cdbpath -run or -delaysetup | |
139 | /// to prevent the parent class from parsing these arguments in DoInit. | |
140 | ||
141 | if (strcmp(arg, "-cdbpath") == 0 or strcmp(arg, "-run") == 0 or | |
142 | strcmp(arg, "-delaysetup") == 0) | |
143 | { | |
144 | return true; | |
145 | } | |
146 | else | |
147 | { | |
148 | return false; | |
149 | } | |
150 | } | |
151 | ||
152 | ||
153 | int AliHLTMUONRawDataHistoComponent::DoInit(int argc, const char** argv) | |
154 | { | |
155 | /// Inherited from AliHLTComponent. | |
156 | /// Parses the command line parameters and initialises the component. | |
157 | ||
158 | HLTInfo("Initialising dHLT raw data histogrammer component."); | |
159 | ||
160 | // Inherit the parents functionality. | |
161 | int result = AliHLTMUONProcessor::DoInit(argc, argv); | |
162 | if (result != 0) return result; | |
163 | ||
164 | fLastPublishTime = fCurrentEventTime = -1; | |
a35ef960 | 165 | fPublishDelay = 0; |
a63da6d6 | 166 | bool pubDelaySet = false; |
7989fd8e | 167 | fSuppressEmptyHists = false; |
168 | fProcessDataEventsOnly = false; | |
a35ef960 | 169 | fClearAfterPublish = false; |
7989fd8e | 170 | fTrackerDecoder.TryRecover(false); |
171 | fTriggerDecoder.TryRecover(false); | |
a63da6d6 | 172 | |
173 | for (int i = 0; i < argc; i++) | |
174 | { | |
175 | if (ArgumentAlreadyHandled(i, argv[i])) continue; | |
176 | ||
177 | if (strcmp(argv[i], "-pubdelay") == 0) | |
178 | { | |
179 | if (pubDelaySet) | |
180 | { | |
181 | HLTWarning("The publishing delay value was already specified." | |
182 | " Will replace previous value given by -pubdelay." | |
183 | ); | |
184 | } | |
185 | ||
186 | if (argc <= i+1) | |
187 | { | |
188 | HLTError("The value for the publishing delay was not specified."); | |
189 | return -EINVAL; | |
190 | } | |
191 | ||
192 | char* cpErr = NULL; | |
193 | double num = strtod(argv[i+1], &cpErr); | |
194 | if (cpErr == NULL or *cpErr != '\0' or num < 0) | |
195 | { | |
196 | HLTError("Cannot convert '%s' to a positive floating point number.", | |
197 | argv[i+1] | |
198 | ); | |
199 | return -EINVAL; | |
200 | } | |
201 | fPublishDelay = num; | |
202 | pubDelaySet = true; | |
203 | ||
204 | i++; | |
205 | continue; | |
206 | } | |
7989fd8e | 207 | |
208 | if (strcmp(argv[i], "-noemptyhists") == 0) | |
209 | { | |
210 | fSuppressEmptyHists = true; | |
211 | continue; | |
212 | } | |
213 | ||
214 | if (strcmp(argv[i], "-onlydataevents") == 0) | |
215 | { | |
216 | fProcessDataEventsOnly = true; | |
217 | continue; | |
218 | } | |
219 | ||
a35ef960 | 220 | if (strcmp(argv[i], "-clearafterpub") == 0) |
221 | { | |
222 | fClearAfterPublish = true; | |
223 | continue; | |
224 | } | |
225 | ||
7989fd8e | 226 | if (strcmp(argv[i], "-tryrecover") == 0) |
227 | { | |
228 | fTrackerDecoder.TryRecover(true); | |
229 | fTriggerDecoder.TryRecover(true); | |
230 | continue; | |
231 | } | |
a63da6d6 | 232 | |
233 | HLTError("Unknown option '%s'.", argv[i]); | |
234 | return -EINVAL; | |
235 | } | |
236 | ||
237 | try | |
238 | { | |
239 | char name[256]; | |
240 | char title[1024]; | |
241 | ||
242 | // Do not add to current directory to prevent memory leak warning. | |
243 | // We will not be leaking any memory if we dont add to the directory. | |
244 | TH1::AddDirectory(kFALSE); | |
245 | ||
246 | for (int i = 0; i < 22; i++) | |
247 | { | |
248 | AliHLTInt32_t equipId = AliHLTMUONUtils::DDLNumberToEquipId(i); | |
249 | sprintf(name, "rawDataErrors_%d", equipId); | |
250 | sprintf(title, "Distribution of errors found in raw data from DDL %d.", equipId); | |
a63da6d6 | 251 | fErrorHist[i] = new TH1D(name, title, 40, 0.5, 40.5); |
252 | fErrorHist[i]->SetXTitle("Error code"); | |
253 | fErrorHist[i]->SetYTitle("Number of errors"); | |
254 | } | |
255 | for (int i = 0; i < 20; i++) | |
256 | { | |
257 | AliHLTInt32_t equipId = AliHLTMUONUtils::DDLNumberToEquipId(i); | |
258 | sprintf(name, "manuDistrib_%d", equipId); | |
259 | sprintf(title, "Distribution of MANUs containing raw data in DDL %d.", equipId); | |
260 | fManuHist[i] = new TH1D(name, title, 2048, -0.5, 2047.5); | |
261 | fManuHist[i]->SetXTitle("MANU number (as seen in raw data)"); | |
262 | fManuHist[i]->SetYTitle("Number of signals read."); | |
263 | sprintf(name, "signalDistrib_%d", equipId); | |
264 | sprintf(title, "Distribution of signals in raw data from DDL %d.", equipId); | |
265 | fSignalHist[i] = new TH1D(name, title, 4096, -0.5, 4095.5); | |
266 | fSignalHist[i]->SetXTitle("Channels"); | |
267 | fSignalHist[i]->SetYTitle("dN/dChannel"); | |
268 | } | |
269 | } | |
270 | catch (const std::bad_alloc&) | |
271 | { | |
272 | HLTError("Could not allocate more memory for histogram objects."); | |
273 | FreeObjects(); | |
274 | return -ENOMEM; | |
275 | } | |
276 | ||
277 | return 0; | |
278 | } | |
279 | ||
280 | ||
281 | int AliHLTMUONRawDataHistoComponent::DoDeinit() | |
282 | { | |
283 | /// Inherited from AliHLTComponent. Performs a cleanup of the component. | |
284 | /// Will delete all histogram objects. | |
285 | ||
286 | HLTInfo("Deinitialising dHLT raw data histogrammer component."); | |
287 | fCurrentEventTime = -1; | |
288 | FreeObjects(); | |
289 | return 0; | |
290 | } | |
291 | ||
292 | ||
293 | int AliHLTMUONRawDataHistoComponent::DoEvent( | |
1d8ae082 | 294 | const AliHLTComponentEventData& evtData, |
295 | AliHLTComponentTriggerData& trigData | |
a63da6d6 | 296 | ) |
297 | { | |
298 | /// Inherited from AliHLTProcessor. | |
299 | /// Processes the new event data and generates summary histograms. | |
300 | ||
7989fd8e | 301 | if (fProcessDataEventsOnly and not IsDataEvent()) return 0; // Only process data events. |
302 | ||
a63da6d6 | 303 | fCurrentEventTime = TTimeStamp().AsDouble(); |
304 | ||
305 | const AliHLTComponentBlockData* block = GetFirstInputBlock(AliHLTMUONConstants::DDLRawDataType()); | |
306 | for ( ; block != NULL; block = GetNextInputBlock()) | |
307 | { | |
308 | HLTDebug("Handling block with fDataType = '%s', fPtr = %p," | |
309 | " fSize = %u bytes and fSpecification = 0x%8.8X.", | |
310 | DataType2Text(block->fDataType).c_str(), block->fPtr, | |
311 | block->fSize, block->fSpecification | |
312 | ); | |
313 | ||
314 | if (AliHLTMUONUtils::IsTrackerDDL(block->fSpecification)) | |
315 | { | |
1d8ae082 | 316 | bool decodeOk = ProcessTrackerDDL(block); |
317 | if (not decodeOk and DumpDataOnError()) | |
318 | { | |
319 | DumpEvent(evtData, trigData); | |
320 | return -EIO; | |
321 | } | |
a63da6d6 | 322 | } |
323 | else if (AliHLTMUONUtils::IsTriggerDDL(block->fSpecification)) | |
324 | { | |
1d8ae082 | 325 | bool decodeOk = ProcessTriggerDDL(block); |
326 | if (not decodeOk and DumpDataOnError()) | |
327 | { | |
328 | DumpEvent(evtData, trigData); | |
329 | return -EIO; | |
330 | } | |
a63da6d6 | 331 | } |
332 | else | |
333 | { | |
334 | HLTError("Received a raw data block with an invalid specification of" | |
335 | " 0x%8.8X. Expected raw data only from one DDL and not multiple" | |
336 | " DDLs as indicated by the specification.", | |
337 | block->fSpecification | |
338 | ); | |
339 | } | |
340 | } | |
341 | ||
7989fd8e | 342 | // See if 'fPublishDelay' number of seconds has elapsed or this is the first event, |
343 | // in that case publish the histograms. Do not publish histograms that are empty | |
344 | // if the fSuppressEmptyHists flag is set. | |
345 | if (fLastPublishTime == -1 or fCurrentEventTime - fLastPublishTime >= fPublishDelay) | |
a63da6d6 | 346 | { |
347 | for (int i = 0; i < 22; i++) | |
348 | { | |
7989fd8e | 349 | if (fSuppressEmptyHists and fErrorHist[i]->GetEntries() == 0) continue; |
a63da6d6 | 350 | PushBack(fErrorHist[i], |
351 | AliHLTMUONConstants::HistogramDataType(), | |
352 | AliHLTMUONUtils::DDLNumberToSpec(i) | |
353 | ); | |
a35ef960 | 354 | // If requested, clear histogram when published. |
355 | if (fClearAfterPublish) fErrorHist[i]->Reset("M"); | |
a63da6d6 | 356 | } |
357 | for (int i = 0; i < 20; i++) | |
358 | { | |
359 | AliHLTUInt32_t spec = AliHLTMUONUtils::DDLNumberToSpec(i); | |
7989fd8e | 360 | if (not (fSuppressEmptyHists and fManuHist[i]->GetEntries() == 0)) |
361 | { | |
362 | PushBack(fManuHist[i], AliHLTMUONConstants::HistogramDataType(), spec); | |
a35ef960 | 363 | if (fClearAfterPublish) fManuHist[i]->Reset("M"); |
7989fd8e | 364 | } |
365 | if (not (fSuppressEmptyHists and fSignalHist[i]->GetEntries() == 0)) | |
366 | { | |
367 | PushBack(fSignalHist[i], AliHLTMUONConstants::HistogramDataType(), spec); | |
a35ef960 | 368 | if (fClearAfterPublish) fSignalHist[i]->Reset("M"); |
7989fd8e | 369 | } |
a63da6d6 | 370 | } |
371 | fLastPublishTime = fCurrentEventTime; | |
372 | } | |
373 | ||
374 | return 0; | |
375 | } | |
376 | ||
377 | ||
1d8ae082 | 378 | bool AliHLTMUONRawDataHistoComponent::ProcessTrackerDDL(const AliHLTComponentBlockData* block) |
a63da6d6 | 379 | { |
380 | /// Processes a raw data block from the tracker stations. | |
381 | ||
382 | AliHLTInt32_t ddl = AliHLTMUONUtils::SpecToDDLNumber(block->fSpecification); | |
383 | assert(0 <= ddl and ddl < 20); | |
384 | ||
7989fd8e | 385 | fTrackerDecoder.GetHandler().ErrorHist(fErrorHist[ddl]); |
386 | fTrackerDecoder.GetHandler().ManuHist(fManuHist[ddl]); | |
387 | fTrackerDecoder.GetHandler().SignalHist(fSignalHist[ddl]); | |
388 | ||
a63da6d6 | 389 | if (block->fSize >= sizeof(AliRawDataHeader)) |
390 | { | |
7989fd8e | 391 | AliHLTUInt8_t* payload = reinterpret_cast<AliHLTUInt8_t*>(block->fPtr) |
392 | + sizeof(AliRawDataHeader); | |
393 | UInt_t payloadSize = UInt_t(block->fSize) - sizeof(AliRawDataHeader); | |
1d8ae082 | 394 | return fTrackerDecoder.Decode(payload, payloadSize); |
a63da6d6 | 395 | } |
396 | else | |
397 | { | |
398 | HLTError("Received a raw data block that is too short to be valid." | |
399 | " Its size is only %d bytes", | |
400 | block->fSize | |
401 | ); | |
402 | fErrorHist[ddl]->Fill(40); | |
403 | } | |
1d8ae082 | 404 | return false; |
a63da6d6 | 405 | } |
406 | ||
407 | ||
1d8ae082 | 408 | bool AliHLTMUONRawDataHistoComponent::ProcessTriggerDDL(const AliHLTComponentBlockData* block) |
a63da6d6 | 409 | { |
410 | /// Processes a raw data block from the trigger stations. | |
411 | ||
a63da6d6 | 412 | AliHLTInt32_t ddl = AliHLTMUONUtils::SpecToDDLNumber(block->fSpecification); |
413 | assert(21 <= ddl and ddl < 22); | |
414 | ||
7989fd8e | 415 | fTriggerDecoder.GetHandler().ErrorHist(fErrorHist[ddl]); |
416 | ||
a63da6d6 | 417 | if (block->fSize >= sizeof(AliRawDataHeader)) |
418 | { | |
419 | AliRawDataHeader* header = reinterpret_cast<AliRawDataHeader*>(block->fPtr); | |
7989fd8e | 420 | AliHLTUInt8_t* payload = reinterpret_cast<AliHLTUInt8_t*>(header+1); |
421 | UInt_t payloadSize = UInt_t(block->fSize) - sizeof(AliRawDataHeader); | |
422 | bool scalarEvent = ((header->GetL1TriggerMessage() & 0x1) == 0x1); | |
1d8ae082 | 423 | return fTriggerDecoder.Decode(payload, payloadSize, scalarEvent); |
a63da6d6 | 424 | } |
425 | else | |
426 | { | |
427 | HLTError("Received a raw data block that is too short to be valid." | |
428 | " Its size is only %d bytes", | |
429 | block->fSize | |
430 | ); | |
431 | fErrorHist[ddl]->Fill(40); | |
432 | } | |
1d8ae082 | 433 | return false; |
a63da6d6 | 434 | } |
435 | ||
436 | ||
437 | void AliHLTMUONRawDataHistoComponent::FreeObjects() | |
438 | { | |
439 | /// Deletes all the histogram objects that were allocated. | |
440 | ||
441 | for (int i = 0; i < 22; i++) | |
442 | { | |
443 | if (fErrorHist[i] != NULL) | |
444 | { | |
445 | delete fErrorHist[i]; | |
446 | fErrorHist[i] = NULL; | |
447 | } | |
448 | } | |
449 | for (int i = 0; i < 20; i++) | |
450 | { | |
451 | if (fManuHist[i] != NULL) | |
452 | { | |
453 | delete fManuHist[i]; | |
454 | fManuHist[i] = NULL; | |
455 | } | |
456 | if (fSignalHist[i] != NULL) | |
457 | { | |
458 | delete fSignalHist[i]; | |
459 | fSignalHist[i] = NULL; | |
460 | } | |
461 | } | |
462 | } | |
7989fd8e | 463 |