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