]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/RCU/AliHLTAltroChannelSelectorComponent.cxx
Ped and Gain DAs update.
[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
21     @date   
22     @brief  A filter/selective readout component for Altro data.
23 */
24
25 // see header file for class documentation
26 // or
27 // refer to README to build package
28 // or
29 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
30
31 #include <cassert>
32 #include "AliHLTAltroChannelSelectorComponent.h"
33 #include "AliAltroDecoder.h"
34 #include "AliAltroData.h"
35 #include "AliAltroBunch.h"
36 #include "TMath.h"
37
38 /** ROOT macro for the implementation of ROOT specific class methods */
39 ClassImp(AliHLTAltroChannelSelectorComponent)
40
41 AliHLTAltroChannelSelectorComponent::AliHLTAltroChannelSelectorComponent()
42   :
43   AliHLTProcessor(),
44   fSkipCorrupted(true),
45   fTalkative(false),
46   fStartTimeBin(0),
47   fEndTimeBin(0),
48   fSignalThreshold(0),
49   fRMSThreshold(0)
50 {
51   // see header file for class documentation
52   // or
53   // refer to README to build package
54   // or
55   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
56 }
57
58 AliHLTAltroChannelSelectorComponent::~AliHLTAltroChannelSelectorComponent()
59 {
60   // see header file for class documentation
61 }
62
63 const char* AliHLTAltroChannelSelectorComponent::GetComponentID()
64 {
65   // see header file for class documentation
66   return "AltroChannelSelector";
67 }
68
69 void AliHLTAltroChannelSelectorComponent::GetInputDataTypes(AliHLTComponentDataTypeList& list)
70 {
71   // see header file for class documentation
72   list.clear();
73   list.push_back(kAliHLTDataTypeDDLRaw|kAliHLTDataOriginTPC);
74   list.push_back(kAliHLTDataTypeHwAddr16);
75 }
76
77 AliHLTComponentDataType AliHLTAltroChannelSelectorComponent::GetOutputDataType()
78 {
79   // see header file for class documentation
80   return kAliHLTDataTypeDDLRaw|kAliHLTDataOriginTPC;
81 }
82
83 void AliHLTAltroChannelSelectorComponent::GetOutputDataSize(unsigned long& constBase, double& inputMultiplier)
84 {
85   // see header file for class documentation
86   constBase=0;
87   inputMultiplier=1.0;
88 }
89
90 AliHLTComponent* AliHLTAltroChannelSelectorComponent::Spawn()
91 {
92   // see header file for class documentation
93   return new AliHLTAltroChannelSelectorComponent;
94 }
95
96 int AliHLTAltroChannelSelectorComponent::DoInit(int argc, const char** argv)
97 {
98   // see header file for class documentation
99   int iResult=0;
100   TString argument="";
101   bool bMissingParam=0;
102   char* cpErr=NULL;
103   int i=0;
104   for (; i<argc && iResult>=0; i++) {
105     cpErr=NULL;
106     argument=argv[i];
107     if (argument.IsNull()) continue;
108
109     // -skip-corrupted, just for backward compatibility, not announced
110     if (argument.CompareTo("-skip-corrupted")==0) {
111       fSkipCorrupted=true;
112
113     // -keep-corrupted
114     } else if (argument.CompareTo("-keep-corrupted")==0) {
115       fSkipCorrupted=false;
116
117     // -talkative
118     } else if (argument.CompareTo("-talkative")==0) {
119       fTalkative=true;
120
121     // -start-timebin
122     } else if (argument.CompareTo("-start-timebin")==0) {
123       if ((bMissingParam=(++i>=argc))) break;
124       fStartTimeBin = strtoul( argv[i], &cpErr ,0);
125       if ( *cpErr ) break;
126
127     // -end-timebin
128     } else if (argument.CompareTo("-end-timebin")==0) {
129       if ((bMissingParam=(++i>=argc))) break;
130       fEndTimeBin = strtoul( argv[i], &cpErr ,0);
131       if ( *cpErr ) break;
132
133     // -signal-threshold
134     } else if (argument.CompareTo("-signal-threshold")==0) {
135       if ((bMissingParam=(++i>=argc))) break;
136       fSignalThreshold = strtoul( argv[i], &cpErr ,0);
137       if ( *cpErr ) break;
138
139     // -rms-threshold
140     } else if (argument.CompareTo("-rms-threshold")==0) {
141       if ((bMissingParam=(++i>=argc))) break;
142       fRMSThreshold = strtoul( argv[i], &cpErr ,0);
143       if ( *cpErr ) break;
144
145     } else {
146       HLTError("unknown argument %s", argument.Data());
147       iResult=-EINVAL;
148     }
149   }
150
151   if (cpErr && *cpErr) {
152     HLTError("Cannot convert specifier '%s' for argument '%s'", argv[i], argument.Data());
153     iResult=-EINVAL;
154   } else if (bMissingParam) {
155     HLTError("missing parameter for argument %s", argument.Data());
156     iResult=-EINVAL;
157   }
158
159   return iResult;
160 }
161
162 int AliHLTAltroChannelSelectorComponent::DoDeinit()
163 {
164   // see header file for class documentation
165   return 0;
166 }
167
168 int AliHLTAltroChannelSelectorComponent::DoEvent(const AliHLTComponentEventData& evtData,
169                                                  const AliHLTComponentBlockData* blocks, 
170                                                  AliHLTComponentTriggerData& /*trigData*/,
171                                                  AliHLTUInt8_t* outputPtr, 
172                                                  AliHLTUInt32_t& size,
173                                                  AliHLTComponentBlockDataList& outputBlocks )
174 {
175   // see header file for class documentation
176   int iResult=0;
177   const int cdhSize=32;
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
188   AliAltroDecoder* decoder=NULL;
189   for (pDesc=GetFirstInputBlock(kAliHLTDataTypeDDLRaw); pDesc!=NULL; pDesc=GetNextInputBlock(), blockno++) {
190     iResult=0;
191     if (pDesc->fSize<=32) {
192       continue;
193     }
194
195     // search for the active pad information
196     AliHLTUInt16_t* pActiveHwAddressArray=NULL;
197     int iArraySize=0;
198     if (fSignalThreshold==0 && fRMSThreshold==0) {
199     for (int i=0; i<(int)evtData.fBlockCnt; i++ ) {
200       // search for selection data of hw address type
201       // which matches the data specification of the block
202       if (blocks[i].fDataType==kAliHLTDataTypeHwAddr16 && blocks[i].fSpecification==pDesc->fSpecification) {
203         pActiveHwAddressArray=reinterpret_cast<AliHLTUInt16_t*>(blocks[i].fPtr);
204         iArraySize=blocks[i].fSize/sizeof(AliHLTUInt16_t);
205         break;
206       }
207     }
208     if (pActiveHwAddressArray==NULL) {
209       HLTWarning("no block of type %s for specification 0x%08x available, data block skipped", 
210                  DataType2Text(kAliHLTDataTypeHwAddr16).c_str(), 
211                  pDesc->fSpecification);
212       break;
213     }
214     }
215
216     if (decoder) delete decoder;
217     decoder=new AliAltroDecoder;
218     if (decoder->SetMemory(reinterpret_cast<UChar_t*>(pDesc->fPtr), pDesc->fSize)<0) {
219       HLTWarning("corrupted data block: initialization of decoder failed for block: %s specification %#x size %d",
220                  DataType2Text(pDesc->fDataType).c_str(), pDesc->fSpecification, pDesc->fSize);
221       iResult=-EFAULT;
222     } else {
223       if (decoder->Decode()) {
224         HLTDebug("init decoder %p size %d", pDesc->fPtr,pDesc->fSize);
225       } else {
226         HLTWarning("corrupted data block: decoding failed for raw data block: %s specification %#x size %d",
227                    DataType2Text(pDesc->fDataType).c_str(), pDesc->fSpecification, pDesc->fSize);
228         iResult=-EFAULT;
229       }
230     }
231
232     unsigned int rcuTrailerLength=0;
233     if (iResult>=0 &&
234         ((rcuTrailerLength=decoder->GetRCUTrailerSize())==0 ||
235          rcuTrailerLength>pDesc->fSize-cdhSize)) {
236       if (rcuTrailerLength>0) {
237         HLTWarning("corrupted data block: RCU trailer length exceeds buffer size");
238       } else {
239         HLTWarning("corrupted data block: RCU trailer of zero length"); 
240       }
241       iResult=-EFAULT;
242     }
243
244     if (iResult<0) {
245       // TODO: here the trigger has to come into play. It is up to
246       // policy if a corrupted data block should be kept (original
247       // DDL) or discarded. In any case, the block is not going to
248       // be part of HLTOUT
249       HLTWarning("skipping corrupted data block for event %lu, data block %s 0x%08x", evtData.fEventID,
250                  DataType2Text(pDesc->fDataType).c_str(), pDesc->fSpecification);
251       iResult=0;
252       continue;
253     }
254
255     int iSelected=0;
256     int iTotal=0;
257     int iCorrupted=0;
258     AliHLTUInt32_t iOutputSize=0;
259     AliHLTUInt32_t iNofAltro40=0;
260     AliHLTUInt32_t iCapacity=size;
261     AliAltroData channel;
262     AliAltroBunch altrobunch;
263
264     // first add the RCU trailer
265     AliHLTUInt8_t* pSrc=reinterpret_cast<AliHLTUInt8_t*>(pDesc->fPtr);
266     pSrc+=pDesc->fSize-rcuTrailerLength;
267     if ((iResult=CopyBlockToEnd(outputPtr, iCapacity, iOutputSize, pSrc, rcuTrailerLength))>=0) {
268       assert(iResult==(int)rcuTrailerLength);
269       iOutputSize+=rcuTrailerLength;
270     } else {
271       HLTError("failed to write RCU trailer of length %d for block %d, too little space in output buffer?", rcuTrailerLength, blockno);
272       iResult=-ENOSPC;
273       break;
274     }
275
276     while (decoder->NextChannel(&channel) && iResult>=0) {
277       iTotal++;
278
279       int hwAddress=channel.GetHadd();
280       if (fSignalThreshold!=0) {
281         // treshold by adc counts
282         unsigned int sumSignals=0;
283         unsigned int maxSignal=0;
284         unsigned int nofSignals=0;
285         while(channel.NextBunch(&altrobunch)){
286           const UInt_t *bunchData=altrobunch.GetData();
287           unsigned int time=altrobunch.GetStartTimeBin();
288           for(Int_t i=0;i<altrobunch.GetBunchSize();i++){
289             if(bunchData[i]>0){// disregarding 0 data.
290               if(time+i>=fStartTimeBin && time+i<=fEndTimeBin){
291                 sumSignals+=bunchData[i];
292                 if (maxSignal<bunchData[i]) maxSignal=bunchData[i];
293                 nofSignals++;
294               }
295             }
296           }
297         }
298         if (nofSignals==0 || maxSignal<=(sumSignals/nofSignals)+fSignalThreshold) {
299           continue;
300         }
301
302       } else if (fRMSThreshold!=0) {
303         // treshold by adc counts
304         unsigned int sumSignals=0;
305         unsigned int maxSignal=0;
306         unsigned int nofSignals=0;
307         while(channel.NextBunch(&altrobunch)){
308           const UInt_t *bunchData=altrobunch.GetData();
309           unsigned int time=altrobunch.GetStartTimeBin();
310           for(Int_t i=0;i<altrobunch.GetBunchSize();i++){
311             if(bunchData[i]>0){// disregarding 0 data.
312               if(time+i>=fStartTimeBin && time+i<=fEndTimeBin){
313                 sumSignals+=bunchData[i]*bunchData[i];
314                 if (maxSignal<bunchData[i]) maxSignal=bunchData[i];
315                 nofSignals++;
316               }
317             }
318           }
319         }
320         if (nofSignals==0 || maxSignal<=TMath::Sqrt(sumSignals/nofSignals)*fRMSThreshold) {
321           continue;
322         }
323         
324       } else {
325       int active=0;
326       for (active=0; active<iArraySize; active++) {
327         if (pActiveHwAddressArray[active]==(AliHLTUInt16_t)hwAddress) {
328           break;
329         }
330       }
331       if (active>=iArraySize) {
332         HLTDebug("ALTRO block %#x (%d) discarded (inactive)", hwAddress, hwAddress);
333         continue;
334       }
335       }
336
337       // no of 10 bit words is without the fill words to fill complete 40 bit words
338       // in addition, align to complete 40 bit words (the '+3')
339       // also, the 5 bytes of the Altro trailer must be added to get the full size
340       int channelSize=((channel.GetDataSize()+3)/4)*5;
341       if (channelSize==0) {
342         if (fTalkative) HLTWarning("skipping zero length channel (hw address %d)", hwAddress);
343         iCorrupted++;
344         continue;
345       }
346       channelSize+=5;
347       HLTDebug("ALTRO block hwAddress 0x%08x (%d) selected (active), size %d", hwAddress, hwAddress, channelSize);
348
349       if ((iResult=decoder->CopyBackward(outputPtr, iCapacity-iOutputSize))>=0) {
350         if (channelSize == iResult) {
351           if (channelSize%5 == 0) {
352             iNofAltro40+=channelSize/5;
353             iOutputSize+=channelSize;
354           } else {
355             if (fTalkative) HLTWarning("corrupted ALTRO channel: incomplete 40 bit word (channel hw address %d)", hwAddress);
356             iCorrupted++;
357             continue;
358           }
359         } else {
360           if (fTalkative) HLTWarning("internal error: failed to copy full channel: %d out of %d bytes (hw address %d)", iResult, channelSize, hwAddress);
361           iCorrupted++;
362           continue;
363         }
364       } else {
365         if (fTalkative) HLTError("failed to write ALTRO channel of length %d for block %d  (hw address %d)", channelSize, blockno, hwAddress);
366         // corrupted channel, but keep going
367         iCorrupted++;
368         iResult=0;
369         continue;
370       }
371       iSelected++;
372     }
373     if (iResult>=0) {
374       // write the common data header
375       if ((iResult=CopyBlockToEnd(outputPtr, iCapacity, iOutputSize, pDesc->fPtr, cdhSize))>=0) {
376         assert(iResult==cdhSize);
377         iOutputSize+=cdhSize;
378
379         // set new length of the data block
380         AliHLTUInt32_t* pCdhSize=reinterpret_cast<AliHLTUInt32_t*>(outputPtr+iCapacity-iOutputSize);
381         *pCdhSize=iOutputSize;
382
383         // set number of Altro words
384         AliHLTUInt32_t* pNofAltro40=reinterpret_cast<AliHLTUInt32_t*>(outputPtr+iCapacity-rcuTrailerLength);
385         *pNofAltro40=iNofAltro40;
386
387         // insert block descriptor
388         AliHLTComponentBlockData bd;
389         FillBlockData(bd);
390         bd.fOffset=iCapacity-iOutputSize;
391         bd.fSize=iOutputSize;
392         bd.fDataType=pDesc->fDataType;
393         bd.fSpecification=pDesc->fSpecification;
394         outputBlocks.push_back(bd);
395         iCapacity-=iOutputSize;
396       } else {
397         HLTError("failed to write CDH of length %d for block %d", cdhSize, blockno);
398         break;
399       }
400     }
401     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);
402   }
403   if (decoder) delete decoder;
404
405   if (iResult<0) {
406     outputBlocks.clear();
407   }
408
409   // all data blocks need to be moved to the beginning of the
410   // buffer because PubSub is not able to handle data blocks entirely
411   // at the end of the buffer. The problem is that the component always
412   // indicates to use the full size of the buffer
413   if (outputBlocks.size()>0) {
414     int offset=outputBlocks.back().fOffset;
415     size-=offset;
416     memmove(outputPtr, outputPtr+offset, size);
417     for (AliHLTComponentBlockDataList::iterator block=outputBlocks.begin();
418          block!=outputBlocks.end();
419          block++) {
420       block->fOffset-=offset;
421     }
422   }
423
424   return iResult;
425 }
426
427 int AliHLTAltroChannelSelectorComponent::CopyBlockToEnd(AliHLTUInt8_t* pTgt, unsigned capacity, unsigned position, void* pSrc, unsigned size)
428 {
429   int iResult=0;
430   if (pTgt==NULL || pSrc==NULL) return -EINVAL;
431   if (capacity-position<size) return -ENOSPC;
432   
433   memcpy(pTgt+(capacity-position-size), pSrc, size);
434   iResult=size;
435   
436   return iResult;
437 }