1 /**************************************************************************
2 * This file is property of and copyright by the ALICE HLT Project *
3 * All rights reserved. *
6 * Artur Szostak <artursz@iafrica.com> *
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 **************************************************************************/
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.
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.
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"
35 #include "AliHLTCDHWrapper.h"
36 #include "TTimeStamp.h"
44 // Helper type for memory allocation.
45 typedef const AliHLTMUONMansoTrackStruct* AliHLTMUONMansoTrackStructP;
48 ClassImp(AliHLTMUONRawDataHistoComponent);
51 AliHLTMUONRawDataHistoComponent::AliHLTMUONRawDataHistoComponent() :
52 AliHLTMUONProcessor(),
56 fCurrentEventTime(-1),
58 fSuppressEmptyHists(false),
59 fProcessDataEventsOnly(false),
60 fClearAfterPublish(false)
62 /// Default constructor initialises all histogram object pointers to NULL.
64 for (int i = 0; i < 22; i++)
68 for (int i = 0; i < 20; i++)
71 fSignalHist[i] = NULL;
74 fTrackerDecoder.ExitOnError(false);
75 fTrackerDecoder.TryRecover(false);
76 fTrackerDecoder.SendDataOnParityError(true);
77 fTrackerDecoder.AutoDetectTrailer(true);
78 fTrackerDecoder.CheckForTrailer(true);
80 fTriggerDecoder.ExitOnError(false);
81 fTriggerDecoder.TryRecover(false);
82 fTriggerDecoder.AutoDetectScalars(false);
86 AliHLTMUONRawDataHistoComponent::~AliHLTMUONRawDataHistoComponent()
88 /// Default destructor deletes all histogram objects if they are still allocated.
93 const char* AliHLTMUONRawDataHistoComponent::GetComponentID()
95 /// Inherited from AliHLTComponent. Returns the component ID.
97 return AliHLTMUONConstants::RawDataHistogrammerId();
101 void AliHLTMUONRawDataHistoComponent::GetInputDataTypes(AliHLTComponentDataTypeList& list)
103 /// Inherited from AliHLTProcessor. Returns the list of expected input data types.
105 assert( list.empty() );
106 list.push_back( AliHLTMUONConstants::DDLRawDataType() );
110 AliHLTComponentDataType AliHLTMUONRawDataHistoComponent::GetOutputDataType()
112 /// Inherited from AliHLTComponent. Returns kAliHLTHistogramDataTypeID.
114 return AliHLTMUONConstants::HistogramDataType();
118 void AliHLTMUONRawDataHistoComponent::GetOutputDataSize(
119 unsigned long& constBase, double& inputMultiplier
122 /// Inherited from AliHLTComponent. Returns an estimate of the expected output data size.
124 constBase = (sizeof(TH1D)+50*sizeof(double))*22 + (sizeof(TH1D)+1024*4*sizeof(double))*20*2;
129 AliHLTComponent* AliHLTMUONRawDataHistoComponent::Spawn()
131 /// Inherited from AliHLTComponent. Creates a new object instance.
133 return new AliHLTMUONRawDataHistoComponent;
137 bool AliHLTMUONRawDataHistoComponent::IgnoreArgument(const char* arg) const
139 /// Return true if the argument is one of -cdbpath -run or -delaysetup
140 /// to prevent the parent class from parsing these arguments in DoInit.
142 if (strcmp(arg, "-cdbpath") == 0 or strcmp(arg, "-run") == 0 or
143 strcmp(arg, "-delaysetup") == 0)
154 int AliHLTMUONRawDataHistoComponent::DoInit(int argc, const char** argv)
156 /// Inherited from AliHLTComponent.
157 /// Parses the command line parameters and initialises the component.
159 HLTInfo("Initialising dHLT raw data histogrammer component.");
161 // Inherit the parents functionality.
162 int result = AliHLTMUONProcessor::DoInit(argc, argv);
163 if (result != 0) return result;
165 fLastPublishTime = fCurrentEventTime = -1;
167 bool pubDelaySet = false;
168 fSuppressEmptyHists = false;
169 fProcessDataEventsOnly = false;
170 fClearAfterPublish = false;
171 fTrackerDecoder.TryRecover(false);
172 fTriggerDecoder.TryRecover(false);
174 for (int i = 0; i < argc; i++)
176 if (ArgumentAlreadyHandled(i, argv[i])) continue;
178 if (strcmp(argv[i], "-pubdelay") == 0)
182 HLTWarning("The publishing delay value was already specified."
183 " Will replace previous value given by -pubdelay."
189 HLTError("The value for the publishing delay was not specified.");
194 double num = strtod(argv[i+1], &cpErr);
195 if (cpErr == NULL or *cpErr != '\0' or num < 0)
197 HLTError("Cannot convert '%s' to a positive floating point number.",
209 if (strcmp(argv[i], "-noemptyhists") == 0)
211 fSuppressEmptyHists = true;
215 if (strcmp(argv[i], "-onlydataevents") == 0)
217 fProcessDataEventsOnly = true;
221 if (strcmp(argv[i], "-clearafterpub") == 0)
223 fClearAfterPublish = true;
227 if (strcmp(argv[i], "-tryrecover") == 0)
229 fTrackerDecoder.TryRecover(true);
230 fTriggerDecoder.TryRecover(true);
234 HLTError("Unknown option '%s'.", argv[i]);
243 // Do not add to current directory to prevent memory leak warning.
244 // We will not be leaking any memory if we dont add to the directory.
245 TH1::AddDirectory(kFALSE);
247 for (int i = 0; i < 22; i++)
249 AliHLTInt32_t equipId = AliHLTMUONUtils::DDLNumberToEquipId(i);
250 sprintf(name, "rawDataErrors_%d", equipId);
251 sprintf(title, "Distribution of errors found in raw data from DDL %d.", equipId);
252 fErrorHist[i] = new TH1D(name, title, 40, 0.5, 40.5);
253 fErrorHist[i]->SetXTitle("Error code");
254 fErrorHist[i]->SetYTitle("Number of errors");
256 for (int i = 0; i < 20; i++)
258 AliHLTInt32_t equipId = AliHLTMUONUtils::DDLNumberToEquipId(i);
259 sprintf(name, "manuDistrib_%d", equipId);
260 sprintf(title, "Distribution of MANUs containing raw data in DDL %d.", equipId);
261 fManuHist[i] = new TH1D(name, title, 2048, -0.5, 2047.5);
262 fManuHist[i]->SetXTitle("MANU number (as seen in raw data)");
263 fManuHist[i]->SetYTitle("Number of signals read.");
264 sprintf(name, "signalDistrib_%d", equipId);
265 sprintf(title, "Distribution of signals in raw data from DDL %d.", equipId);
266 fSignalHist[i] = new TH1D(name, title, 4096, -0.5, 4095.5);
267 fSignalHist[i]->SetXTitle("Channels");
268 fSignalHist[i]->SetYTitle("dN/dChannel");
271 catch (const std::bad_alloc&)
273 HLTError("Could not allocate more memory for histogram objects.");
282 int AliHLTMUONRawDataHistoComponent::DoDeinit()
284 /// Inherited from AliHLTComponent. Performs a cleanup of the component.
285 /// Will delete all histogram objects.
287 HLTInfo("Deinitialising dHLT raw data histogrammer component.");
288 fCurrentEventTime = -1;
294 int AliHLTMUONRawDataHistoComponent::DoEvent(
295 const AliHLTComponentEventData& evtData,
296 AliHLTComponentTriggerData& trigData
299 /// Inherited from AliHLTProcessor.
300 /// Processes the new event data and generates summary histograms.
302 if (fProcessDataEventsOnly and not IsDataEvent()) return 0; // Only process data events.
304 fCurrentEventTime = TTimeStamp().AsDouble();
306 const AliHLTComponentBlockData* block = GetFirstInputBlock(AliHLTMUONConstants::DDLRawDataType());
307 for ( ; block != NULL; block = GetNextInputBlock())
309 HLTDebug("Handling block with fDataType = '%s', fPtr = %p,"
310 " fSize = %u bytes and fSpecification = 0x%8.8X.",
311 DataType2Text(block->fDataType).c_str(), block->fPtr,
312 block->fSize, block->fSpecification
315 if (AliHLTMUONUtils::IsTrackerDDL(block->fSpecification))
317 bool decodeOk = ProcessTrackerDDL(block);
318 if (not decodeOk and DumpDataOnError())
320 DumpEvent(evtData, trigData);
324 else if (AliHLTMUONUtils::IsTriggerDDL(block->fSpecification))
326 bool decodeOk = ProcessTriggerDDL(block);
327 if (not decodeOk and DumpDataOnError())
329 DumpEvent(evtData, trigData);
335 HLTError("Received a raw data block with an invalid specification of"
336 " 0x%8.8X. Expected raw data only from one DDL and not multiple"
337 " DDLs as indicated by the specification.",
338 block->fSpecification
343 // See if 'fPublishDelay' number of seconds has elapsed or this is the first event,
344 // in that case publish the histograms. Do not publish histograms that are empty
345 // if the fSuppressEmptyHists flag is set.
346 if (fLastPublishTime == -1 or fCurrentEventTime - fLastPublishTime >= fPublishDelay)
348 for (int i = 0; i < 22; i++)
350 if (fSuppressEmptyHists and fErrorHist[i]->GetEntries() == 0) continue;
351 PushBack(fErrorHist[i],
352 AliHLTMUONConstants::HistogramDataType(),
353 AliHLTMUONUtils::DDLNumberToSpec(i)
355 // If requested, clear histogram when published.
356 if (fClearAfterPublish) fErrorHist[i]->Reset("M");
358 for (int i = 0; i < 20; i++)
360 AliHLTUInt32_t spec = AliHLTMUONUtils::DDLNumberToSpec(i);
361 if (not (fSuppressEmptyHists and fManuHist[i]->GetEntries() == 0))
363 PushBack(fManuHist[i], AliHLTMUONConstants::HistogramDataType(), spec);
364 if (fClearAfterPublish) fManuHist[i]->Reset("M");
366 if (not (fSuppressEmptyHists and fSignalHist[i]->GetEntries() == 0))
368 PushBack(fSignalHist[i], AliHLTMUONConstants::HistogramDataType(), spec);
369 if (fClearAfterPublish) fSignalHist[i]->Reset("M");
372 fLastPublishTime = fCurrentEventTime;
379 bool AliHLTMUONRawDataHistoComponent::ProcessTrackerDDL(const AliHLTComponentBlockData* block)
381 /// Processes a raw data block from the tracker stations.
383 AliHLTInt32_t ddl = AliHLTMUONUtils::SpecToDDLNumber(block->fSpecification);
384 assert(0 <= ddl and ddl < 20);
386 fTrackerDecoder.GetHandler().ErrorHist(fErrorHist[ddl]);
387 fTrackerDecoder.GetHandler().ManuHist(fManuHist[ddl]);
388 fTrackerDecoder.GetHandler().SignalHist(fSignalHist[ddl]);
390 AliHLTCDHWrapper cdh(block->fPtr);
391 if (block->fSize >= sizeof(AliRawDataHeader) &&
392 block->fSize >= cdh.GetHeaderSize()) // in case if cdh v3
394 AliHLTUInt8_t* payload = reinterpret_cast<AliHLTUInt8_t*>(block->fPtr)
395 + cdh.GetHeaderSize();
396 UInt_t payloadSize = UInt_t(block->fSize) - cdh.GetHeaderSize();
397 return fTrackerDecoder.Decode(payload, payloadSize);
401 HLTError("Received a raw data block that is too short to be valid."
402 " Its size is only %d bytes",
405 fErrorHist[ddl]->Fill(40);
411 bool AliHLTMUONRawDataHistoComponent::ProcessTriggerDDL(const AliHLTComponentBlockData* block)
413 /// Processes a raw data block from the trigger stations.
415 AliHLTInt32_t ddl = AliHLTMUONUtils::SpecToDDLNumber(block->fSpecification);
416 assert(21 <= ddl and ddl < 22);
418 fTriggerDecoder.GetHandler().ErrorHist(fErrorHist[ddl]);
420 AliHLTCDHWrapper cdh(block->fPtr);
421 if (block->fSize >= sizeof(AliRawDataHeader) &&
422 block->fSize >= cdh.GetHeaderSize()) // in case if cdh v3)
424 AliHLTUInt8_t* payload = reinterpret_cast<AliHLTUInt8_t*>(block->fPtr);
425 payload += cdh.GetHeaderSize();
426 UInt_t payloadSize = UInt_t(block->fSize) - cdh.GetHeaderSize();
427 bool scalarEvent = ((cdh.GetL1TriggerMessage() & 0x1) == 0x1);
428 return fTriggerDecoder.Decode(payload, payloadSize, scalarEvent);
432 HLTError("Received a raw data block that is too short to be valid."
433 " Its size is only %d bytes",
436 fErrorHist[ddl]->Fill(40);
442 void AliHLTMUONRawDataHistoComponent::FreeObjects()
444 /// Deletes all the histogram objects that were allocated.
446 for (int i = 0; i < 22; i++)
448 if (fErrorHist[i] != NULL)
450 delete fErrorHist[i];
451 fErrorHist[i] = NULL;
454 for (int i = 0; i < 20; i++)
456 if (fManuHist[i] != NULL)
461 if (fSignalHist[i] != NULL)
463 delete fSignalHist[i];
464 fSignalHist[i] = NULL;