1 #ifndef AliHLTMUONRAWDATAHISTOCOMPONENT_H
2 #define AliHLTMUONRAWDATAHISTOCOMPONENT_H
3 /* This file is property of and copyright by the ALICE HLT Project *
4 * ALICE Experiment at CERN, All rights reserved. *
5 * See cxx source for full Copyright notice */
10 /// @file AliHLTMUONRawDataHistoComponent.h
11 /// @author Artur Szostak <artursz@iafrica.com>
12 /// @date 30 April 2008
13 /// @brief Declaration of a component to generate basic monitoring histograms of raw data.
16 #include "AliHLTMUONProcessor.h"
17 #include "AliHLTMUONDataTypes.h"
18 #include "AliMUONTrackerDDLDecoder.h"
19 #include "AliMUONTriggerDDLDecoder.h"
23 * @class AliHLTMUONRawDataHistoComponent
24 * @brief Dimuon HLT component for generating basic monitoring histograms for raw data.
26 * This component is useful for performing basic monitoring tasks on the raw data
27 * from the muon spectrometer. It will try and decode the data and histogram the
28 * following information:
29 * \li The distribution of signals per DDL.
30 * \li The number of ADC values found per MANU for each DDL.
31 * \li The error codes found by the decoders while trying to decode the data for each DDL.
33 * <h2>General properties:</h2>
35 * Component ID: \b MUONRawDataHistogrammer <br>
36 * Library: \b libAliHLTMUON.so <br>
37 * Input Data Types: AliHLTMUONConstants::DDLRawDataType() = "DDL_RAW :MUON" <br>
38 * Output Data Types: AliHLTMUONConstants::HistogramDataType() = "ROOTHIST:MUON" <br>
40 * <h2>Mandatory arguments:</h2>
43 * <h2>Optional arguments:</h2>
44 * \li -pubdelay <i>delay</i> <br>
45 * Indicates the number of seconds to wait between publishing the histograms.
46 * The default is zero seconds. <i>delay</i> must be a positive floating point
48 * \li -noemptyhists <br>
49 * If indicated then any histograms that are empty will not be published.
50 * By default all events are processed. <br>
51 * \li -onlydataevents <br>
52 * If indicated then only data events are processed.
53 * By default all events are processed. <br>
54 * \li -clearafterpub <br>
55 * If specified then all the internal histograms are cleared after they are
56 * published, so they will not accumilate statistics over the whole run.
57 * This is off by default. <br>
58 * \li -tryrecover <br>
59 * This is a special option to the raw data decoder to turn on logic which will
60 * try and recover from corrupt raw DDL data. This is off by default. <br>
62 * <h2>Standard configuration:</h2>
63 * There is no special configuration for this component.
65 * <h2>Default CDB entries:</h2>
68 * <h2>Performance:</h2>
69 * A few milliseconds per event.
71 * <h2>Memory consumption:</h2>
72 * Minimal, under 1 MBytes.
74 * <h2>Output size:</h2>
77 * @ingroup alihlt_dimuon_component
79 class AliHLTMUONRawDataHistoComponent : public AliHLTMUONProcessor
82 AliHLTMUONRawDataHistoComponent();
83 virtual ~AliHLTMUONRawDataHistoComponent();
85 // Public functions to implement the AliHLTProcessor interface.
86 // These functions are required for the registration process.
87 virtual const char* GetComponentID();
88 virtual void GetInputDataTypes(AliHLTComponentDataTypeList& list);
89 virtual AliHLTComponentDataType GetOutputDataType();
90 virtual void GetOutputDataSize(unsigned long& constBase, double& inputMultiplier);
91 virtual AliHLTComponent* Spawn();
95 // Protected functions to implement the AliHLTProcessor interface.
96 // These functions provide initialization as well as the actual processing
97 // capabilities of the component.
98 virtual int DoInit(int argc, const char** argv);
99 virtual bool IgnoreArgument(const char* arg) const;
100 virtual int DoDeinit();
102 const AliHLTComponentEventData& evtData,
103 AliHLTComponentTriggerData& trigData
106 using AliHLTProcessor::DoEvent;
110 class AliDecoderHandler : public AliHLTLogging
114 AliDecoderHandler() : AliHLTLogging(), fErrorHist(NULL) {}
115 virtual ~AliDecoderHandler() {}
117 /// Returns the error codes histogram.
118 TH1D* ErrorHist() const { return fErrorHist; }
120 /// Sets the error codes histogram.
121 void ErrorHist(TH1D* hist) { fErrorHist = hist; }
125 /// Fills the error histogram with the given error code.
126 void FillErrorHist(Int_t code);
128 TH1D* fErrorHist; /// Histogram of error codes found.
132 // Do not allow copying of this object.
134 AliDecoderHandler(const AliDecoderHandler& obj);
136 AliDecoderHandler& operator = (const AliDecoderHandler& obj);
139 class AliTrackerDecoderHandler :
140 public AliMUONTrackerDDLDecoderEventHandler, public AliDecoderHandler
143 AliTrackerDecoderHandler() :
144 AliMUONTrackerDDLDecoderEventHandler(),
150 virtual ~AliTrackerDecoderHandler() {}
152 /// Returns the signals per MANU histogram.
153 TH1D* ManuHist() const { return fManuHist; }
155 /// Sets the signals per MANU histogram.
156 void ManuHist(TH1D* hist) { fManuHist = hist; }
158 /// Returns the signals histogram.
159 TH1D* SignalHist() const { return fSignalHist; }
161 /// Sets the signals histogram.
162 void SignalHist(TH1D* hist) { fSignalHist = hist; }
164 // Methods inherited from AliMUONTrackerDDLDecoderEventHandler:
166 /// Called for each new data word found.
167 void OnData(UInt_t data, bool /*parityError*/);
169 /// Fills the fErrorHist histogram with the error code received.
170 void OnError(ErrorCode code, const void* /*location*/) { FillErrorHist(Int_t(code)); }
174 // Do not allow copying of this object.
176 AliTrackerDecoderHandler(const AliTrackerDecoderHandler& obj);
178 AliTrackerDecoderHandler& operator = (const AliTrackerDecoderHandler& obj);
180 TH1D* fManuHist; /// Histogram of signal distribution per MANU.
181 TH1D* fSignalHist; /// Histogram of the ADC signal distribution.
184 class AliTriggerDecoderHandler :
185 public AliMUONTriggerDDLDecoderEventHandler, public AliDecoderHandler
188 AliTriggerDecoderHandler() :
189 AliMUONTriggerDDLDecoderEventHandler(),
193 virtual ~AliTriggerDecoderHandler() {}
195 // Methods inherited from AliMUONTriggerDDLDecoderEventHandler:
197 /// Fills the fErrorHist histogram with the error code received.
198 void OnError(ErrorCode code, const void* /*location*/) { FillErrorHist(Int_t(code)); }
202 // Do not allow copying of this object.
204 AliTriggerDecoderHandler(const AliTriggerDecoderHandler& obj);
206 AliTriggerDecoderHandler& operator = (const AliTriggerDecoderHandler& obj);
209 // Do not allow copying of this class.
210 AliHLTMUONRawDataHistoComponent(const AliHLTMUONRawDataHistoComponent& /*obj*/);
211 AliHLTMUONRawDataHistoComponent& operator = (const AliHLTMUONRawDataHistoComponent& /*obj*/);
214 * Decodes the tracker DDL data block and fills the histograms.
215 * \param block The data block to decode.
216 * \returns true if the block could be decoded and false if there was an error in the data.
218 bool ProcessTrackerDDL(const AliHLTComponentBlockData* block);
221 * Decodes the trigger DDL data block and fills the histograms.
222 * \param block The data block to decode.
223 * \returns true if the block could be decoded and false if there was an error in the data.
225 bool ProcessTriggerDDL(const AliHLTComponentBlockData* block);
228 * Deletes all the histograms and resets the pointers.
232 AliMUONTrackerDDLDecoder<AliTrackerDecoderHandler> fTrackerDecoder; // Raw data decoder for the tracker data.
233 AliMUONTriggerDDLDecoder<AliTriggerDecoderHandler> fTriggerDecoder; // Raw data decoder for the trigger data.
235 double fLastPublishTime; /// Timestamp for the last time we published data (seconds).
236 double fCurrentEventTime; /// Timestamp for the current event being processed (seconds).
237 double fPublishDelay; /// Delay in second to wait between publishing data.
238 TH1D* fErrorHist[22]; /// Histograms for error codes per DDL.
239 TH1D* fManuHist[20]; /// Histograms for MANU distributions per DDL.
240 TH1D* fSignalHist[20]; /// Histograms for signal distributions per DDL.
241 bool fSuppressEmptyHists; /// Flag indicating if empty histograms should be published or not.
242 bool fProcessDataEventsOnly; /// Flag indicating if only data events should be processed.
243 bool fClearAfterPublish; /// Inficates if the histograms should be reset after being published.
245 ClassDef(AliHLTMUONRawDataHistoComponent, 0); // Trigger decision component for the dimuon HLT.
248 //-----------------------------------------------------------------------------
250 inline void AliHLTMUONRawDataHistoComponent::AliDecoderHandler::FillErrorHist(Int_t code)
252 /// Fills the error code into the error code histogram.
254 assert(fErrorHist != NULL);
255 Int_t mincode = Int_t( fErrorHist->GetXaxis()->GetBinCenter(1) );
256 Int_t maxcode = Int_t( fErrorHist->GetXaxis()->GetBinCenter(fErrorHist->GetNbinsX()) );
257 if (code < mincode or maxcode < code)
259 HLTError("Filling an error code which is out of range."
260 " Received code %d, but expected it to be in the range [%d..%d]",
261 int(code), mincode, maxcode
264 fErrorHist->Fill(code);
268 inline void AliHLTMUONRawDataHistoComponent::AliTrackerDecoderHandler::OnData(
269 UInt_t data, bool /*parityError*/
272 /// Fills the signals histogram and also the signal per MANU histogram.
274 UInt_t minadc = UInt_t( fSignalHist->GetXaxis()->GetBinCenter(1) );
275 UInt_t maxadc = UInt_t( fSignalHist->GetXaxis()->GetBinCenter(fSignalHist->GetNbinsX()) );
276 UInt_t minmanu = UInt_t( fManuHist->GetXaxis()->GetBinCenter(1) );
277 UInt_t maxmanu = UInt_t( fManuHist->GetXaxis()->GetBinCenter(fManuHist->GetNbinsX()) );
279 UShort_t manuId; UChar_t channelId; UShort_t adc;
280 UnpackADC(data, manuId, channelId, adc);
282 if (adc < minadc or maxadc < adc)
284 HLTError("Filling a signal value which is out of range. Received ADC value"
285 " of %d channels, but expected it to be in the range [%d..%d]",
286 int(adc), minadc, maxadc
289 fSignalHist->Fill(adc);
291 if (manuId < minmanu or maxmanu < manuId)
293 HLTError("Filling a MANU ID value which is out of range. Received"
294 " value of %d, but expected it to be in the range [%d..%d]",
295 int(manuId), minmanu, maxmanu
298 fManuHist->Fill(manuId);
301 #endif // AliHLTMUONRAWDATAHISTOCOMPONENT_H