]> git.uio.no Git - u/mrichter/AliRoot.git/blame - HLT/PHOS/AliHLTPHOSRawAnalyzerComponent.cxx
Made the fEquipmentID variable constant
[u/mrichter/AliRoot.git] / HLT / PHOS / AliHLTPHOSRawAnalyzerComponent.cxx
CommitLineData
cbab66dd 1/**************************************************************************
99388135 2 * This file is property of and copyright by the Experimental Nuclear *
3 * Physics Group, Dep. of Physics *
4 * University of Oslo, Norway, 2007 *
5 * *
6 * Author: Per Thomas Hille <perthi@fys.uio.no> for the ALICE HLT Project.*
cbab66dd 7 * Contributors are mentioned in the code where appropriate. *
99388135 8 * Please report bugs to perthi@fys.uio.no *
cbab66dd 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
2b7cf4fb 19
20/// @class AliHLTPHOSRawAnalyzerComponent
21/// Base class of PHOS HLT online raw analysis component.
22/// The class provides a common interface for the implementation of PHOS
23/// HLT raw data
24/// processors components. The class is intended for processing of
25/// arrays of raw data samples to evaluate energy and timing.
26/// The Energy will be given in entities of ADC leves ranging from 0 to
27/// 1023. Timing will be given in entities of samples periods.
28/// Drived clases must implement the fucntions
29/// - @ref GetComponentID
30/// - @ref Spawn
31
32
33
34#include "AliHLTPHOSRawAnalyzer.h"
cbab66dd 35#include "AliHLTPHOSRawAnalyzerComponent.h"
9dfd64cf 36#include "AliRawReaderMemory.h"
0a211711 37#include "AliCaloRawStream.h"
2ca5f55f 38#include <iostream>
bde48b84 39#include "AliHLTPHOSRcuCellEnergyDataStruct.h"
2b7cf4fb 40#include "AliHLTPHOSRcuChannelDataStruct.h"
cbab66dd 41
2ca5f55f 42using namespace std;
43
44
146c463a 45const AliHLTComponentDataType AliHLTPHOSRawAnalyzerComponent::fgkInputDataTypes[]={kAliHLTVoidDataType,{0,"",""}}; //'zero' terminated array
46int AliHLTPHOSRawAnalyzerComponent::fgEventCount = 0;
ee7849e6 47
4df5dd10 48
2b7cf4fb 49
50AliHLTPHOSRawAnalyzerComponent::AliHLTPHOSRawAnalyzerComponent():AliHLTProcessor(), fAnalyzerPtr(0),
51fEquippmentID(0), fModuleID(0), fRcuX(0), fRcuZ(0),fRcuZOffset(0), fRcuXOffset(0),fPrintInfo(kFALSE),fSendChannelData(kFALSE),fPrintInfoFrequncy(1000),
52fPHOSRawStream(0), fRawMemoryReader(0), fOutPtr(0)
cbab66dd 53{
2ca5f55f 54 /**
55 *Default constructor
56 */
cbab66dd 57}
58
2b7cf4fb 59
cbab66dd 60AliHLTPHOSRawAnalyzerComponent::~AliHLTPHOSRawAnalyzerComponent()
61{
2ca5f55f 62 /**
63 *Default destructor
64 */
0a211711 65 if(fRawMemoryReader != 0)
66 {
67 delete fRawMemoryReader;
68 }
69 if(fPHOSRawStream != 0)
70 {
71 delete fPHOSRawStream;
72 }
cbab66dd 73}
74
75
146c463a 76AliHLTPHOSRawAnalyzerComponent::AliHLTPHOSRawAnalyzerComponent(const AliHLTPHOSRawAnalyzerComponent & ) : AliHLTProcessor(), fAnalyzerPtr(0),
2b7cf4fb 77fEquippmentID(0), fModuleID(0), fRcuX(0), fRcuZ(0),fRcuZOffset(0), fRcuXOffset(0),fPrintInfo(kFALSE),fSendChannelData(kFALSE),fPrintInfoFrequncy(1000),
78fPHOSRawStream(0), fRawMemoryReader(0), fOutPtr(0)
cbab66dd 79{
2ca5f55f 80 /**
81 *Copy Constructor
82 */
cbab66dd 83}
84
85int
86AliHLTPHOSRawAnalyzerComponent::Deinit()
87{
2b7cf4fb 88 //See base class for documentation
857f8ed5 89 cout << "Deinit" << endl;
cbab66dd 90 return 0;
91}
92
93int
94AliHLTPHOSRawAnalyzerComponent::DoDeinit()
95{
2b7cf4fb 96 //See base class for documentation
857f8ed5 97 cout << "DoDeinit" << endl;
9dfd64cf 98 Logging(kHLTLogInfo, "HLT", "PHOS", ",AliHLTPHOSRawAnalyzerComponen DoDeinit");
99
0a211711 100 if(fRawMemoryReader !=0)
101 {
102 delete fRawMemoryReader;
103 }
104
105 if(fPHOSRawStream != 0)
106 {
107 delete fPHOSRawStream;
108 }
cbab66dd 109 return 0;
cbab66dd 110
cbab66dd 111}
112
9dfd64cf 113const char*
114AliHLTPHOSRawAnalyzerComponent::GetComponentID()
115{
2b7cf4fb 116 ///Returns the component ID
9dfd64cf 117 return "AliPhosTestRaw";
118}
ee7849e6 119
1c1b3412 120
cbab66dd 121void
ee7849e6 122AliHLTPHOSRawAnalyzerComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list)
cbab66dd 123{
2b7cf4fb 124 //See Base class for documentation
146c463a 125 const AliHLTComponentDataType* pType=fgkInputDataTypes;
ee7849e6 126 while (pType->fID!=0) {
127 list.push_back(*pType);
128 pType++;
129 }
cbab66dd 130}
131
132AliHLTComponentDataType
133AliHLTPHOSRawAnalyzerComponent::GetOutputDataType()
134{
2b7cf4fb 135 //See Base class for documentation
53740333 136 return AliHLTPHOSDefinitions::gkCellEnergyDataType;
cbab66dd 137}
138
139void
9dfd64cf 140AliHLTPHOSRawAnalyzerComponent::GetOutputDataSize(unsigned long& constBase, double& inputMultiplier )
141
cbab66dd 142{
2b7cf4fb 143 //See Base class for documentation
ef408bb3 144 constBase = 30;
2b7cf4fb 145 // inputMultiplier = 0.1;
146 inputMultiplier = 1;
cbab66dd 147}
148
0a211711 149int AliHLTPHOSRawAnalyzerComponent::DoEvent( const AliHLTComponentEventData& evtData, const AliHLTComponentBlockData* blocks,
150 AliHLTComponentTriggerData& trigData, AliHLTUInt8_t* outputPtr,
151 AliHLTUInt32_t& size, vector<AliHLTComponentBlockData>& outputBlocks )
cbab66dd 152{
2b7cf4fb 153 /// Function that proceesses the raw date to give Energy and TOF for the
154 /// Individual cells/crystals.
155
857f8ed5 156 AliHLTUInt8_t tmpMod = 0;
157 AliHLTUInt8_t tmpZ = 0;
158 AliHLTUInt8_t tmpX = 0;
159 AliHLTUInt8_t tmpGain = 0;
bde48b84 160 Int_t sampleCnt = 0;
2947a32c 161 Int_t processedChannels = 0;
bde48b84 162 UInt_t offset = 0;
163 UInt_t mysize = 0;
164 UInt_t tSize = 0;
2bcb5a06 165 Int_t tmpChannelCnt = 0;
432edd34 166 Int_t tmpStartIndex = 0;
53740333 167 AliHLTUInt8_t* outBPtr;
2b7cf4fb 168 unsigned long first;
169 unsigned long last;
53740333 170 outBPtr = outputPtr;
0a211711 171 const AliHLTComponentBlockData* iter = NULL;
172 unsigned long ndx;
cf434398 173
0a211711 174 for( ndx = 0; ndx < evtData.fBlockCnt; ndx++ )
175 {
0a211711 176 iter = blocks+ndx;
53740333 177 mysize = 0;
178 offset = tSize;
179
0a211711 180 if ( iter->fDataType != AliHLTPHOSDefinitions::gkDDLPackedRawDataType )
181 {
0a211711 182 continue;
183 }
ef408bb3 184
bde48b84 185 fRawMemoryReader->SetMemory( reinterpret_cast<UChar_t*>( iter->fPtr ), iter->fSize );
146c463a 186 fAnalyzerPtr->SetData(fTmpChannelData);
a00e1689 187 fOutPtr = (AliHLTPHOSRcuCellEnergyDataStruct*)outBPtr;
bde48b84 188 mysize += sizeof(AliHLTPHOSRcuCellEnergyDataStruct);
a00e1689 189 fOutPtr->fRcuX = fRcuX;
190 fOutPtr->fRcuZ = fRcuZ;
191 fOutPtr->fModuleID = fModuleID;
2bcb5a06 192 tmpChannelCnt = 0;
2b7cf4fb 193
2947a32c 194 while(fPHOSRawStream->Next())
195 {
196 if (fPHOSRawStream->IsNewHWAddress())
197 {
198 if(processedChannels > 0)
199 {
146c463a 200 fAnalyzerPtr->SetData(fTmpChannelData);
201 fAnalyzerPtr->Evaluate(0, sampleCnt);
a00e1689 202 fOutPtr->fValidData[tmpChannelCnt].fGain = tmpGain;
1c1b3412 203 fOutPtr->fValidData[tmpChannelCnt].fZ = tmpZ;
204 fOutPtr->fValidData[tmpChannelCnt].fX = tmpX;
2b7cf4fb 205 fOutPtr->fValidData[tmpChannelCnt].fEnergy = (float)fAnalyzerPtr->GetEnergy();
206 fOutPtr->fValidData[tmpChannelCnt].fTime = (float)fAnalyzerPtr->GetTiming();
432edd34 207 ResetDataPtr(tmpStartIndex, sampleCnt);
2b7cf4fb 208 tmpChannelCnt ++;
68d9caee 209 sampleCnt = 0;
2947a32c 210 }
211
1c1b3412 212 tmpMod = (AliHLTUInt8_t)fPHOSRawStream->GetModule() ;
1c1b3412 213 tmpX =(AliHLTUInt8_t)fPHOSRawStream->GetRow() - fRcuXOffset;
214 tmpZ =(AliHLTUInt8_t)fPHOSRawStream->GetColumn() - fRcuZOffset;
2947a32c 215 tmpGain = fPHOSRawStream->IsLowGain();
216 processedChannels ++;
217 }
432edd34 218
219 if(sampleCnt == 0)
220 {
221 tmpStartIndex = fPHOSRawStream->GetTime();
222 }
223
2947a32c 224 fTmpChannelData[fPHOSRawStream->GetTime()] = fPHOSRawStream->GetSignal();
bde48b84 225 sampleCnt ++;
2947a32c 226 }
2b7cf4fb 227
006ee9b2 228 tmpChannelCnt ++;
146c463a 229 fAnalyzerPtr->SetData(fTmpChannelData);
230 fAnalyzerPtr->Evaluate(0, sampleCnt);
006ee9b2 231 fOutPtr->fValidData[tmpChannelCnt].fGain = tmpGain;
232 fOutPtr->fValidData[tmpChannelCnt].fZ = tmpZ;
233 fOutPtr->fValidData[tmpChannelCnt].fX = tmpX;
146c463a 234 fOutPtr->fValidData[tmpChannelCnt].fEnergy = fAnalyzerPtr->GetEnergy();
235 fOutPtr->fValidData[tmpChannelCnt].fTime = fAnalyzerPtr->GetTiming();
857f8ed5 236 ResetDataPtr(tmpStartIndex, sampleCnt);
006ee9b2 237
2b7cf4fb 238 sampleCnt = 0;
006ee9b2 239
a00e1689 240 fOutPtr->fCnt = tmpChannelCnt;
53740333 241 AliHLTComponentBlockData bd;
242 FillBlockData( bd );
53740333 243 bd.fOffset = offset;
244 bd.fSize = mysize;
53740333 245 bd.fDataType = AliHLTPHOSDefinitions::gkCellEnergyDataType;
9ce19a20 246 bd.fSpecification = 0xFFFFFFFF;
53740333 247 outputBlocks.push_back( bd );
53740333 248 tSize += mysize;
249 outBPtr += mysize;
2b7cf4fb 250
cf434398 251 if( tSize > size )
2b7cf4fb 252 {
253 cout <<"kHLTLogFatal, HLT::AliHLTPHOSRawAnalyzerComponent::DoEvent Too much data Data written over allowed buffer. Amount written:" << tSize << " allowed" << size << endl;
254 Logging( kHLTLogFatal, "HLT::AliHLTPHOSRawAnalyzerComponent::DoEvent", "Too much data",
255 "Data written over allowed buffer. Amount written: %lu, allowed amount: %lu."
256 , tSize, size );
257 return EMSGSIZE;
258 }
53740333 259
2b7cf4fb 260 }
261
146c463a 262 fgEventCount++;
2b7cf4fb 263
264 if(fPrintInfo == kTRUE)
265 {
266 if(fgEventCount%fPrintInfoFrequncy == 0)
267 {
268 cout <<"Analyzing event " << fgEventCount << "for Equippment " << fEquippmentID << endl;
269 }
270 }
53740333 271 size = tSize;
2947a32c 272 return 0;
273}//end DoEvent
0a211711 274
275
2b7cf4fb 276
0a211711 277int
278AliHLTPHOSRawAnalyzerComponent::DoInit( int argc, const char** argv )
279{
2b7cf4fb 280 /// See base class for documentation
281 fSendChannelData = kFALSE;
282 fPrintInfo = kFALSE;
2947a32c 283 Reset();
0a211711 284 fRawMemoryReader = new AliRawReaderMemory();
285 fPHOSRawStream = new AliCaloRawStream(fRawMemoryReader,"PHOS");
1c1b3412 286 fPHOSRawStream->SetOldRCUFormat(kFALSE);
2b7cf4fb 287 int iResult=0;
288 TString argument="";
289 Bool_t isSetEquippmentID = kFALSE;
290
291 for(int i=0; i<argc && iResult>=0; i++)
292 {
293 argument=argv[i];
294
295 if (argument.IsNull())
296 {
297 continue;
298 }
299
300 if (argument.CompareTo("-equipmentID") == 0)
301 {
302 cout << "AliHLTPHOSRawAnalyzerComponent:DoInit argument = -equipmentID " <<endl;
303 if(i+1 <= argc)
304 {
eae02153 305 SetEquippmentID((AliHLTUInt16_t)atoi(argv[i+1]));
2b7cf4fb 306 cout << "AliHLTPHOSRawAnalyzerComponent:DoInit setting equippment ID to " << fEquippmentID <<endl;
307 fRawMemoryReader->SetEquipmentID(fEquippmentID);
2b7cf4fb 308 SetCoordinates(fEquippmentID);
309 isSetEquippmentID = kTRUE;
310 }
311 else
312 {
313 iResult= -1;
314 Logging( kHLTLogFatal, "HLT::AliHLTPHOSRcuHistogramProducerComponent::DoInt( int argc, const char** argv )", "Missing argument",
315 "The argument -equippmentID expects a number");
316 return iResult;
317 }
318 }
319
320
321 if (argument.CompareTo("-datatype") == 0)
322 {
323 if(i+1 <= argc)
324 {
325 argument=argv[i+1];
326 if(argument.CompareTo("channeldata") == 0)
327 {
328 cout << "AliHLTPHOSRawAnalyzerComponent::DoIni setting sendChannelData = kTRUE "<< endl;
329 fSendChannelData = kTRUE;
330 }
331 }
332 }
333
334 if (argument.CompareTo("-printinfo") == 0)
335 {
336 if(i+1 <= argc)
337 {
338 argument=argv[i+1];
339 fPrintInfoFrequncy = atoi(argv[i+1]);
340 fPrintInfo = kTRUE;
341 cout << "AliHLTPHOSRawAnalyzerComponent::DoIni setting printinfo = kTRUE, with update frequency every "<< fPrintInfoFrequncy << "th event" <<endl;
342 }
343 else
344 {
345 cout << "WARNING: asking for event info, but no update frequency is specified, otipn is ignored" << endl;
346 }
347 }
348
349 }
350
351
352 if(isSetEquippmentID == kFALSE)
353 {
354 Logging( kHLTLogFatal, "HLT::AliHLTPHOSRcuHistogramProducerComponent::DoInt( int argc, const char** argv )", "Missing argument",
355 "The argument equippmentID is not set: set it with a component argumet like this: -equippmentID <number>");
356 iResult = -2;
357 }
358
0a211711 359 return 0;
360}
9dfd64cf 361
05be0766 362
2ca5f55f 363void
364AliHLTPHOSRawAnalyzerComponent::DumpData(int gain) const
0a211711 365{
2b7cf4fb 366 //Dumping data to std out
432edd34 367 for(int mod = 0; mod < N_MODULES; mod ++)
2947a32c 368 {
369 printf("\n *********** MODULE %d ************\n", mod);
432edd34 370 for(int row = 0; row < N_ROWS_MOD; row ++)
2947a32c 371 {
432edd34 372 for(int col = 0; col < N_COLUMNS_MOD; col ++)
2947a32c 373 {
374 if( fMaxValues[mod][row][col][0] != 0)
375 {
432edd34 376 cout << fMaxValues[mod][row][col][gain] << "\t";
2947a32c 377 }
378 }
379 }
380 }
381}
9dfd64cf 382
432edd34 383
2ca5f55f 384void
385AliHLTPHOSRawAnalyzerComponent::DumpChannelData(Double_t *data) const
68d9caee 386{
2ca5f55f 387 //shutting up the code checker
388 cout << endl;
432edd34 389 for(int i=0; i< ALTRO_MAX_SAMPLES; i++)
68d9caee 390 {
391 if (data[i] != 0)
392 {
393 cout <<i <<"\t";
394 }
395 }
396 cout << endl;
397
432edd34 398 for(int i=0; i< ALTRO_MAX_SAMPLES; i++)
68d9caee 399 {
400 if (data[i] != 0)
401 {
402 cout <<data[i] <<"\t";
403 }
404 }
405
406 cout << endl;
407}
408
409
2947a32c 410void
411AliHLTPHOSRawAnalyzerComponent::Reset()
412{
2ca5f55f 413 //shutting code checker
432edd34 414 for(int mod = 0; mod < N_MODULES; mod ++)
9dfd64cf 415 {
432edd34 416 for(int row = 0; row < N_ROWS_MOD; row ++)
2947a32c 417 {
432edd34 418 for(int col = 0; col < N_COLUMNS_MOD; col ++)
2947a32c 419 {
432edd34 420 for(int gain = 0; gain < N_GAINS; gain ++ )
2947a32c 421 {
422 fMaxValues[mod][row][col][gain] = 0;
423 }
424 }
425 }
426 }
9dfd64cf 427
2b7cf4fb 428 ResetDataPtr(0, ALTRO_MAX_SAMPLES);
432edd34 429
430} // end Reset
431
2947a32c 432void
432edd34 433AliHLTPHOSRawAnalyzerComponent::ResetDataPtr(int startindex, int sampleCnt)
2947a32c 434{
2ca5f55f 435 //shutting up the code checker
432edd34 436 for(int i = startindex ; i< sampleCnt; i++)
2947a32c 437 {
438 fTmpChannelData[i] = 0;
439 }
cbab66dd 440}
ef408bb3 441
ef408bb3 442void
1c1b3412 443AliHLTPHOSRawAnalyzerComponent::SetEquippmentID(AliHLTUInt16_t id)
ef408bb3 444{
eae02153 445 ///Changing the value of the constant fEquippmentID
446 ///by virue of const_cast as it should only be set once
447 ///and then remain constant. It caannot be set in the class constructor
448 ///because it should be set in the DoInit fucntion.
449 AliHLTUInt16_t &ref = const_cast<AliHLTUInt16_t&>(fEquippmentID);
450 ref = id;
ef408bb3 451}
452
05be0766 453
2b7cf4fb 454const AliHLTUInt16_t
2ca5f55f 455AliHLTPHOSRawAnalyzerComponent::GetEquippmentID() const
53740333 456{
2ca5f55f 457 //shutting up the code checker
53740333 458 return fEquippmentID;
459}
460
53740333 461void
1c1b3412 462AliHLTPHOSRawAnalyzerComponent::SetCoordinates(AliHLTUInt16_t equippmentID)
ef408bb3 463{
2ca5f55f 464 //shutting up the code checker
1c1b3412 465 int rcuIndex = (fEquippmentID - 1792)%N_RCUS_PER_MODULE;
1c1b3412 466 fModuleID = (fEquippmentID -1792 -rcuIndex)/N_RCUS_PER_MODULE;
467
53740333 468 if(rcuIndex == 0)
469 {
470 fRcuX = 0;
cf434398 471 fRcuZ = 0;
53740333 472 }
473
474 if(rcuIndex == 1)
475 {
cf434398 476 fRcuX = 0;
477 fRcuZ = 1;
53740333 478 }
479
480 if(rcuIndex == 2)
481 {
cf434398 482 fRcuX = 1;
483 fRcuZ = 0;
53740333 484 }
485
cf434398 486 if(rcuIndex == 3)
53740333 487 {
488 fRcuX = 1;
cf434398 489 fRcuZ = 1;
53740333 490 }
491
1c1b3412 492 fRcuZOffset = N_ZROWS_RCU*fRcuZ;
493 fRcuXOffset = N_XCOLUMNS_RCU*fRcuX;
cf434398 494
1c1b3412 495 cout <<"********InitInfo************"<< endl;
eae02153 496 cout <<"AliHLTPHOSRawAnalyzerComponent::SetCoordinate casted"<< endl;
1c1b3412 497 cout <<"Equpippment ID =\t"<< fEquippmentID <<endl;
006ee9b2 498 cout <<"Module ID =\t"<< (int)fModuleID<<endl;
499 cout <<"RCUX =\t\t" << (int)fRcuX << endl;
500 cout <<"RCUZ =\t\t" << (int)fRcuZ << endl;
501 cout <<"RcuZOffset = \t" << (int)fRcuZOffset << endl;
502 cout <<"RcuXOffset = \t" << (int)fRcuXOffset << endl << endl;
ef408bb3 503}