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