]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/RCU/AliHLTAltroChannelSelectorComponent.cxx
updated version of the ALTRO channel selector component using AliAltroRawStreamV3
[u/mrichter/AliRoot.git] / HLT / RCU / AliHLTAltroChannelSelectorComponent.cxx
1 // $Id$
2
3 //**************************************************************************
4 //* This file is property of and copyright by the ALICE HLT Project        * 
5 //* ALICE Experiment at CERN, All rights reserved.                         *
6 //*                                                                        *
7 //* Primary Authors: Matthias Richter <Matthias.Richter@ift.uib.no>        *
8 //*                  for The ALICE HLT Project.                            *
9 //*                                                                        *
10 //* Permission to use, copy, modify and distribute this software and its   *
11 //* documentation strictly for non-commercial purposes is hereby granted   *
12 //* without fee, provided that the above copyright notice appears in all   *
13 //* copies and that both the copyright notice and this permission notice   *
14 //* appear in the supporting documentation. The authors make no claims     *
15 //* about the suitability of this software for any purpose. It is          *
16 //* provided "as is" without express or implied warranty.                  *
17 //**************************************************************************
18
19 /// @file   AliHLTAltroChannelSelectorComponent.cxx
20 /// @author Matthias Richter edited by Jason Glyndwr Ulery
21 /// @date   
22 /// @brief  A filter/selective readout component for Altro data.
23 ///
24
25 #include <cassert>
26 #include <memory>
27 #include "AliHLTAltroChannelSelectorComponent.h"
28 #include "AliHLTErrorGuard.h"
29 #include "AliHLTDAQ.h"
30 #include "AliRawReaderMemory.h"
31 #include "AliAltroRawStreamV3.h"
32 #include "TMath.h"
33
34 /** ROOT macro for the implementation of ROOT specific class methods */
35 ClassImp(AliHLTAltroChannelSelectorComponent)
36
37 AliHLTAltroChannelSelectorComponent::AliHLTAltroChannelSelectorComponent()
38   :
39   AliHLTProcessor(),
40   fSkipCorrupted(true),
41   fTalkative(false),
42   fStartTimeBin(0),
43   fEndTimeBin(0),
44   fSignalThreshold(0),
45   fRMSThreshold(0)
46 {
47   // see header file for class documentation
48   // or
49   // refer to README to build package
50   // or
51   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
52 }
53
54 AliHLTAltroChannelSelectorComponent::~AliHLTAltroChannelSelectorComponent()
55 {
56   // see header file for class documentation
57 }
58
59 const char* AliHLTAltroChannelSelectorComponent::GetComponentID()
60 {
61   // see header file for class documentation
62   return "AltroChannelSelector";
63 }
64
65 void AliHLTAltroChannelSelectorComponent::GetInputDataTypes(AliHLTComponentDataTypeList& list)
66 {
67   // see header file for class documentation
68   list.clear();
69   list.push_back(kAliHLTDataTypeDDLRaw|kAliHLTDataOriginTPC);
70   list.push_back(kAliHLTDataTypeHwAddr16);
71 }
72
73 AliHLTComponentDataType AliHLTAltroChannelSelectorComponent::GetOutputDataType()
74 {
75   // see header file for class documentation
76   return kAliHLTDataTypeDDLRaw|kAliHLTDataOriginTPC;
77 }
78
79 void AliHLTAltroChannelSelectorComponent::GetOutputDataSize(unsigned long& constBase, double& inputMultiplier)
80 {
81   // see header file for class documentation
82   constBase=0;
83   inputMultiplier=1.0;
84 }
85
86 AliHLTComponent* AliHLTAltroChannelSelectorComponent::Spawn()
87 {
88   // see header file for class documentation
89   return new AliHLTAltroChannelSelectorComponent;
90 }
91
92 int AliHLTAltroChannelSelectorComponent::DoInit(int argc, const char** argv)
93 {
94   // see header file for class documentation
95   int iResult=0;
96   TString argument="";
97   bool bMissingParam=0;
98   char* cpErr=NULL;
99   int i=0;
100   for (; i<argc && iResult>=0; i++) {
101     cpErr=NULL;
102     argument=argv[i];
103     if (argument.IsNull()) continue;
104
105     // -skip-corrupted, just for backward compatibility, not announced
106     if (argument.CompareTo("-skip-corrupted")==0) {
107       fSkipCorrupted=true;
108
109     // -keep-corrupted
110     } else if (argument.CompareTo("-keep-corrupted")==0) {
111       fSkipCorrupted=false;
112
113     // -talkative
114     } else if (argument.CompareTo("-talkative")==0) {
115       fTalkative=true;
116
117     // -start-timebin
118     } else if (argument.CompareTo("-start-timebin")==0) {
119       if ((bMissingParam=(++i>=argc))) break;
120       fStartTimeBin = strtoul( argv[i], &cpErr ,0);
121       if ( *cpErr ) break;
122
123     // -end-timebin
124     } else if (argument.CompareTo("-end-timebin")==0) {
125       if ((bMissingParam=(++i>=argc))) break;
126       fEndTimeBin = strtoul( argv[i], &cpErr ,0);
127       if ( *cpErr ) break;
128
129     // -signal-threshold
130     } else if (argument.CompareTo("-signal-threshold")==0) {
131       if ((bMissingParam=(++i>=argc))) break;
132       fSignalThreshold = strtoul( argv[i], &cpErr ,0);
133       if ( *cpErr ) break;
134
135     // -rms-threshold
136     } else if (argument.CompareTo("-rms-threshold")==0) {
137       if ((bMissingParam=(++i>=argc))) break;
138       fRMSThreshold = strtoul( argv[i], &cpErr ,0);
139       if ( *cpErr ) break;
140
141     } else {
142       HLTError("unknown argument %s", argument.Data());
143       iResult=-EINVAL;
144     }
145   }
146
147   if (cpErr && *cpErr) {
148     HLTError("Cannot convert specifier '%s' for argument '%s'", argv[i], argument.Data());
149     iResult=-EINVAL;
150   } else if (bMissingParam) {
151     HLTError("missing parameter for argument %s", argument.Data());
152     iResult=-EINVAL;
153   }
154
155   return iResult;
156 }
157
158 int AliHLTAltroChannelSelectorComponent::DoDeinit()
159 {
160   // see header file for class documentation
161   return 0;
162 }
163
164 int AliHLTAltroChannelSelectorComponent::DoEvent(const AliHLTComponentEventData& evtData,
165                                                  const AliHLTComponentBlockData* blocks, 
166                                                  AliHLTComponentTriggerData& /*trigData*/,
167                                                  AliHLTUInt8_t* outputPtr, 
168                                                  AliHLTUInt32_t& size,
169                                                  AliHLTComponentBlockDataList& outputBlocks )
170 {
171   // see header file for class documentation
172   int iResult=0;
173   AliHLTUInt32_t iCapacity=size;
174   size=0;
175
176   const int cdhSize=32;//8 32-bit words so 32 bytes
177   const int rcuSize=36;//9 32-bit words
178
179   if (!IsDataEvent()) {
180     size=0;
181     return 0;
182   }
183
184   // process the DLL input
185   int blockno=0;
186   const AliHLTComponentBlockData* pDesc=NULL;
187   std::auto_ptr<AliRawReaderMemory> pRawReader(new AliRawReaderMemory);
188   if (!pRawReader.get()) return -ENOMEM;
189
190   for (pDesc=GetFirstInputBlock(kAliHLTDataTypeDDLRaw); pDesc!=NULL; pDesc=GetNextInputBlock(), blockno++) {
191     iResult=0;
192     if (pDesc->fSize<=32) {
193       continue;
194     }
195
196     // search for the active pad information
197     AliHLTUInt16_t* pActiveHwAddressArray=NULL;
198     int iArraySize=0;
199     if (fSignalThreshold==0 && fRMSThreshold==0) {
200     for (int i=0; i<(int)evtData.fBlockCnt; i++ ) {
201       // search for selection data of hw address type
202       // which matches the data specification of the block
203       if (blocks[i].fDataType==kAliHLTDataTypeHwAddr16 && blocks[i].fSpecification==pDesc->fSpecification) {
204         pActiveHwAddressArray=reinterpret_cast<AliHLTUInt16_t*>(blocks[i].fPtr);
205         iArraySize=blocks[i].fSize/sizeof(AliHLTUInt16_t);
206         break;
207       }
208     }
209     if (pActiveHwAddressArray==NULL) {
210       HLTWarning("no block of type %s for specification 0x%08x available, data block skipped", 
211                  DataType2Text(kAliHLTDataTypeHwAddr16).c_str(), 
212                  pDesc->fSpecification);
213       break;
214     }
215     }
216
217     pRawReader->Reset();
218     int ddlid=AliHLTDAQ::DdlIDFromHLTBlockData(pDesc->fDataType.fOrigin, pDesc->fSpecification);
219     if (ddlid<0) {
220       HLTError("unable to extract DDL Id for data block %s 0x%08x", DataType2Text(pDesc->fDataType).c_str(), pDesc->fSpecification);
221       continue;
222     }
223
224     if (!pRawReader->AddBuffer((UChar_t*)pDesc->fPtr,pDesc->fSize, ddlid)) {
225       ALIHLTERRORGUARD(1, "can not set up AltroDecoder for data block %s 0x%08x,"
226                        " skipping data block and suppressing further messages",
227                        DataType2Text(pDesc->fDataType).c_str(), pDesc->fSpecification);
228       continue;
229     }
230
231     std::auto_ptr<AliAltroRawStreamV3> altroRawStream(new AliAltroRawStreamV3(pRawReader.get()));
232
233     if (!altroRawStream.get()) {
234       iResult=-ENOMEM;
235       break;
236     }
237
238     altroRawStream->Reset();
239     if (!altroRawStream->NextDDL()) {
240       ALIHLTERRORGUARD(1, "internal error, can not read data from AliRawReaderMemory");
241       continue;
242     }
243
244     unsigned int rcuTrailerLength=rcuSize;
245     if(rcuTrailerLength>pDesc->fSize-cdhSize) HLTWarning("corrupted data block: RCU trailer length exceeds buffer size");
246   
247     if (iResult<0) {
248       // TODO: here the trigger has to come into play. It is up to
249       // policy if a corrupted data block should be kept (original
250       // DDL) or discarded. In any case, the block is not going to
251       // be part of HLTOUT
252       HLTWarning("skipping corrupted data block for event %lu, data block %s 0x%08x", evtData.fEventID,
253                  DataType2Text(pDesc->fDataType).c_str(), pDesc->fSpecification);
254       iResult=0;
255       continue;
256     }
257
258     int iSelected=0;
259     int iTotal=0;
260     int iCorrupted=0;
261     AliHLTUInt32_t iOutputSize=0;
262
263     //First add the common data header (CDH)
264     memcpy(outputPtr+size,pDesc->fPtr,cdhSize);
265     iOutputSize+=cdhSize;
266
267     while (altroRawStream->NextChannel() && iResult>=0) {
268       iTotal++;
269       int hwAddress=altroRawStream->GetHWAddress();
270   
271       if (fSignalThreshold!=0) {
272         // treshold by adc counts
273         unsigned int sumSignals=0;
274         unsigned int maxSignal=0;
275         unsigned int nofSignals=0;
276         while(altroRawStream->NextBunch()){
277           const UShort_t *bunchData=altroRawStream->GetSignals();
278           unsigned int time=altroRawStream->GetStartTimeBin();
279           for(Int_t i=0;i<altroRawStream->GetBunchLength();i++){
280             if(bunchData[i]>0){// disregarding 0 data.
281               if(time+i>=fStartTimeBin && time+i<=fEndTimeBin){
282                 sumSignals+=bunchData[i];
283                 if (maxSignal<bunchData[i]) maxSignal=bunchData[i];
284                 nofSignals++;
285               }
286             }
287           }
288         }
289         if (nofSignals==0 || maxSignal<=(sumSignals/nofSignals)+fSignalThreshold) {
290           continue;
291         }
292
293       } else if (fRMSThreshold!=0) {
294         // treshold by adc counts
295         unsigned int sumSignals=0;
296         unsigned int maxSignal=0;
297         unsigned int nofSignals=0;
298         while(altroRawStream->NextBunch()){
299           const UShort_t *bunchData=altroRawStream->GetSignals();
300           unsigned int time=altroRawStream->GetStartTimeBin();
301           for(Int_t i=0;i<altroRawStream->GetBunchLength();i++){
302             if(bunchData[i]>0){// disregarding 0 data.
303               if(time+i>=fStartTimeBin && time+i<=fEndTimeBin){
304                 sumSignals+=bunchData[i]*bunchData[i];
305                 if (maxSignal<bunchData[i]) maxSignal=bunchData[i];
306                 nofSignals++;
307               }
308             }
309           }
310         }
311         if (nofSignals==0 || maxSignal<=TMath::Sqrt(sumSignals/nofSignals)*fRMSThreshold) {
312           continue;
313         }
314         
315       } else {
316       int active=0;
317       for (active=0; active<iArraySize; active++) {
318         if (pActiveHwAddressArray[active]==(AliHLTUInt16_t)hwAddress) {
319           break;
320         }
321       }
322       if (active>=iArraySize) {
323         HLTDebug("ALTRO block %#x (%d) discarded (inactive)", hwAddress, hwAddress);
324         continue;
325       }
326       }
327
328       // no of 10 bit words is without the fill words to fill complete 40 bit words
329       // in addition, align to complete 40 bit words (the '+3')
330       // also, the 4 bytes of the Altro trailer must be added to get the full size
331       int channelSize=((altroRawStream->GetChannelPayloadSize()+2)/3)*4;
332       if (channelSize==0) {
333         if (fTalkative) HLTWarning("skipping zero length channel (hw address %d)", hwAddress);
334         iCorrupted++;
335         continue;
336       }
337       channelSize+=4;
338       HLTDebug("ALTRO block hwAddress 0x%08x (%d) selected (active), size %d", hwAddress, hwAddress, channelSize);
339
340       if (iCapacity>=(channelSize+iOutputSize)) {
341         const UChar_t *ChannelData=altroRawStream->GetChannelPayload();
342         iResult=channelSize;
343         memcpy(outputPtr+size+iOutputSize,ChannelData,channelSize);
344         if (channelSize == iResult){
345           if (channelSize%4 == 0) {
346             iOutputSize+=channelSize;
347           } else {
348             if (fTalkative) HLTWarning("corrupted ALTRO channel: incomplete 32 bit word (channel hw address %d)", hwAddress);
349             iCorrupted++;
350             continue;
351           }
352         } else {
353           if (fTalkative) HLTWarning("internal error: failed to copy full channel: %d out of %d bytes (hw address %d)", iResult, channelSize, hwAddress);
354           iCorrupted++;
355           continue;
356         }
357       } else {
358         if (fTalkative) HLTError("failed to write ALTRO channel of length %d for block %d  (hw address %d)", channelSize, blockno, hwAddress);
359         // corrupted channel, but keep going
360         iCorrupted++;
361         iResult=0;
362         continue;
363       }  
364       iSelected++;
365     }
366     
367     if (iResult>=0) {
368       //Write the RCU Trailer
369       AliHLTUInt8_t* pSrc=reinterpret_cast<AliHLTUInt8_t*>(pDesc->fPtr);
370       pSrc+=pDesc->fSize-rcuTrailerLength;
371       memcpy(outputPtr+size+iOutputSize,pSrc,rcuTrailerLength);
372       iOutputSize+=rcuTrailerLength;
373
374       //Set new payload length of the new data bloock 
375       AliHLTUInt32_t* pSize=reinterpret_cast<AliHLTUInt32_t*>(outputPtr+size+iOutputSize-rcuTrailerLength);
376       (*pSize)&=~0x3FFFFFF;
377       (*pSize)|=(iOutputSize-rcuTrailerLength-cdhSize+3)/4;//# of 32 bit words in payload
378       
379       //Write DDL length     
380       AliHLTUInt32_t* pCdhSize=reinterpret_cast<AliHLTUInt32_t*>(outputPtr+size);
381       *pCdhSize=iOutputSize;
382
383       // insert block descriptor
384       AliHLTComponentBlockData bd;
385       FillBlockData(bd);
386       bd.fOffset=size;
387       bd.fSize=iOutputSize;
388       bd.fDataType=pDesc->fDataType;
389       bd.fSpecification=pDesc->fSpecification;
390       outputBlocks.push_back(bd);
391       size+=iOutputSize;
392     }
393     if (fTalkative) HLTImportant("data block %d (0x%08x): selected %d out of %d ALTRO channel(s), %d corrupted channels skipped", blockno, pDesc->fSpecification, iSelected, iTotal, iCorrupted);
394   }
395
396   if (iResult<0) {
397     outputBlocks.clear();
398   }
399   return iResult;
400 }