]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/BASE/AliHLTComponent.cxx
Moving GUID calculation to base class and setting unique ID of global decisions to...
[u/mrichter/AliRoot.git] / HLT / BASE / AliHLTComponent.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 //*                  Timm Steinbeck <timm@kip.uni-heidelberg.de>           *
9 //*                  for The ALICE HLT Project.                            *
10 //*                                                                        *
11 //* Permission to use, copy, modify and distribute this software and its   *
12 //* documentation strictly for non-commercial purposes is hereby granted   *
13 //* without fee, provided that the above copyright notice appears in all   *
14 //* copies and that both the copyright notice and this permission notice   *
15 //* appear in the supporting documentation. The authors make no claims     *
16 //* about the suitability of this software for any purpose. It is          *
17 //* provided "as is" without express or implied warranty.                  *
18 //**************************************************************************
19
20 //  @file   AliHLTComponent.cxx
21 //  @author Matthias Richter, Timm Steinbeck
22 //  @date   
23 //  @brief  Base class implementation for HLT components. */
24 //  @note   The class is both used in Online (PubSub) and Offline (AliRoot)
25 //          context
26
27
28 #if __GNUC__>= 3
29 using namespace std;
30 #endif
31
32 //#include "AliHLTStdIncludes.h"
33 #include "AliHLTComponent.h"
34 #include "AliHLTComponentHandler.h"
35 #include "AliHLTMessage.h"
36 #include "AliHLTCTPData.h"
37 #include "AliRawDataHeader.h"
38 #include "TString.h"
39 #include "TMath.h"
40 #include "TObjArray.h"
41 #include "TObjectTable.h"
42 #include "TClass.h"
43 #include "TStopwatch.h"
44 #include "TFormula.h"
45 #include "TUUID.h"
46 #include "TMD5.h"
47 #include "TRandom3.h"
48 #include "AliHLTMemoryFile.h"
49 #include "AliHLTMisc.h"
50 #include <cassert>
51 #include <ctime>
52 #include <stdint.h>
53
54 /**
55  * default compression level for ROOT objects
56  */
57 #define ALIHLTCOMPONENT_DEFAULT_OBJECT_COMPRESSION 5
58 #define ALIHLTCOMPONENT_STATTIME_SCALER 1000000
59
60 /** ROOT macro for the implementation of ROOT specific class methods */
61 ClassImp(AliHLTComponent);
62
63 /** stopwatch macro using the stopwatch guard */
64 #define ALIHLTCOMPONENT_STOPWATCH(type) AliHLTStopwatchGuard swguard(fpStopwatches!=NULL?reinterpret_cast<TStopwatch*>(fpStopwatches->At((int)type)):NULL)
65 //#define ALIHLTCOMPONENT_STOPWATCH(type) 
66
67 /** stopwatch macro for operations of the base class */
68 #define ALIHLTCOMPONENT_BASE_STOPWATCH() ALIHLTCOMPONENT_STOPWATCH(kSWBase)
69 /** stopwatch macro for operations of the detector algorithm (DA) */
70 #define ALIHLTCOMPONENT_DA_STOPWATCH() ALIHLTCOMPONENT_STOPWATCH(kSWDA)
71
72 AliHLTComponent::AliHLTComponent()
73   :
74   fEnvironment(),
75   fCurrentEvent(0),
76   fEventCount(-1),
77   fFailedEvents(0),
78   fCurrentEventData(),
79   fpInputBlocks(NULL),
80   fCurrentInputBlock(-1),
81   fSearchDataType(kAliHLTVoidDataType),
82   fClassName(),
83   fpInputObjects(NULL),
84   fpOutputBuffer(NULL),
85   fOutputBufferSize(0),
86   fOutputBufferFilled(0),
87   fOutputBlocks(),
88   fpStopwatches(new TObjArray(kSWTypeCount)),
89   fMemFiles(),
90   fpRunDesc(NULL),
91   fCDBSetRunNoFunc(false),
92   fChainId(),
93   fChainIdCrc(0),
94   fpBenchmark(NULL),
95   fFlags(0),
96   fEventType(gkAliEventTypeUnknown),
97   fComponentArgs(),
98   fEventDoneData(NULL),
99   fEventDoneDataSize(0),
100   fCompressionLevel(ALIHLTCOMPONENT_DEFAULT_OBJECT_COMPRESSION)
101   , fLastObjectSize(0)
102   , fpCTPData(NULL)
103   , fPushbackPeriod(0)
104   , fLastPushBackTime(-1)
105 {
106   // see header file for class documentation
107   // or
108   // refer to README to build package
109   // or
110   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
111   memset(&fEnvironment, 0, sizeof(AliHLTAnalysisEnvironment));
112   if (fgpComponentHandler)
113     fgpComponentHandler->ScheduleRegister(this);
114   //SetLocalLoggingLevel(kHLTLogDefault);
115 }
116
117 AliHLTComponent::~AliHLTComponent()
118 {
119   // see header file for function documentation
120   if (fpBenchmark) delete fpBenchmark;
121   fpBenchmark=NULL;
122
123   CleanupInputObjects();
124   if (fpStopwatches!=NULL) delete fpStopwatches;
125   fpStopwatches=NULL;
126   AliHLTMemoryFilePList::iterator element=fMemFiles.begin();
127   while (element!=fMemFiles.end()) {
128     if (*element) {
129       if ((*element)->IsClosed()==0) {
130         HLTWarning("memory file has not been closed, possible data loss or incomplete buffer");
131         // close but do not flush as we dont know whether the buffer is still valid
132         (*element)->CloseMemoryFile(0);
133       }
134       delete *element;
135       *element=NULL;
136     }
137     element++;
138   }
139   if (fpRunDesc) {
140     delete fpRunDesc;
141     fpRunDesc=NULL;
142   }
143   if (fEventDoneData)
144     delete [] reinterpret_cast<AliHLTUInt8_t*>( fEventDoneData );
145   fEventDoneData=NULL;
146
147   if (fpCTPData) {
148     delete fpCTPData;
149   }
150   fpCTPData=NULL;
151 }
152
153 AliHLTComponentHandler* AliHLTComponent::fgpComponentHandler=NULL;
154
155 int AliHLTComponent::SetGlobalComponentHandler(AliHLTComponentHandler* pCH, int bOverwrite) 
156 {
157   // see header file for function documentation
158   int iResult=0;
159   if (fgpComponentHandler==NULL || bOverwrite!=0)
160     fgpComponentHandler=pCH;
161   else
162     iResult=-EPERM;
163   return iResult;
164 }
165
166 int AliHLTComponent::UnsetGlobalComponentHandler() 
167 {
168   // see header file for function documentation
169   return SetGlobalComponentHandler(NULL,1);
170 }
171
172 int AliHLTComponent::SetComponentEnvironment(const AliHLTAnalysisEnvironment* comenv, void* environParam)
173 {
174   // see header file for function documentation
175   HLTLogKeyword(fChainId.c_str());
176   int iResult=0;
177   if (comenv) {
178     memset(&fEnvironment, 0, sizeof(AliHLTAnalysisEnvironment));
179     memcpy(&fEnvironment, comenv, comenv->fStructSize<sizeof(AliHLTAnalysisEnvironment)?comenv->fStructSize:sizeof(AliHLTAnalysisEnvironment));
180     fEnvironment.fStructSize=sizeof(AliHLTAnalysisEnvironment);
181     fEnvironment.fParam=environParam;
182   }
183   return iResult;
184 }
185
186 int AliHLTComponent::Init(const AliHLTAnalysisEnvironment* comenv, void* environParam, int argc, const char** argv )
187 {
188   // see header file for function documentation
189   HLTLogKeyword(fChainId.c_str());
190   int iResult=0;
191   if (comenv) {
192     SetComponentEnvironment(comenv, environParam);
193   }
194   fPushbackPeriod=0;
195   fLastPushBackTime=-1;
196
197   fComponentArgs="";
198   const char** pArguments=NULL;
199   int iNofChildArgs=0;
200   TString argument="";
201   int bMissingParam=0;
202   if (argc>0) {
203     pArguments=new const char*[argc];
204     if (pArguments) {
205       for (int i=0; i<argc && iResult>=0; i++) {
206         if (fComponentArgs.size()>0) fComponentArgs+=" ";
207         fComponentArgs+=argv[i];
208         argument=argv[i];
209         if (argument.IsNull()) continue;
210
211         // benchmark
212         if (argument.CompareTo("-benchmark")==0) {
213
214           // -loglevel=
215         } else if (argument.BeginsWith("-loglevel=")) {
216           TString parameter=argument.ReplaceAll("-loglevel=", "");
217           parameter.Remove(TString::kLeading, ' '); // remove all blanks
218           if (parameter.BeginsWith("0x") &&
219               parameter.Replace(0,2,"",0).IsHex()) {
220             unsigned int loglevel=kHLTLogNone;
221             sscanf(parameter.Data(),"%x", &loglevel);
222             SetLocalLoggingLevel((AliHLTComponentLogSeverity)loglevel);
223           } else {
224             HLTError("wrong parameter for argument %s, hex number expected", argument.Data());
225             iResult=-EINVAL;
226           }
227           // -object-compression=
228         } else if (argument.BeginsWith("-object-compression=")) {
229           argument.ReplaceAll("-object-compression=", "");
230           if (argument.IsDigit()) {
231             fCompressionLevel=argument.Atoi();
232             if (fCompressionLevel<0 || fCompressionLevel>9) {
233               HLTWarning("invalid compression level %d, setting to default %d", fCompressionLevel, ALIHLTCOMPONENT_DEFAULT_OBJECT_COMPRESSION);
234               fCompressionLevel=ALIHLTCOMPONENT_DEFAULT_OBJECT_COMPRESSION;
235             }
236           } else {
237             HLTError("wrong parameter for argument -object-compression, number expected");
238           }
239           // -pushback-period=
240         } else if (argument.BeginsWith("-pushback-period=")) {
241           argument.ReplaceAll("-pushback-period=", "");
242           if (argument.IsDigit()) {
243             fPushbackPeriod=argument.Atoi();
244           } else {
245             HLTError("wrong parameter for argument -pushback-period, number expected");
246           }
247           // -disable-component-stat
248         } else if (argument.CompareTo("-disable-component-stat")==0) {
249           fFlags|=kDisableComponentStat;
250         } else {
251           pArguments[iNofChildArgs++]=argv[i];
252         }
253       }
254     } else {
255       iResult=-ENOMEM;
256     }
257   }
258   if (bMissingParam) {
259     HLTError("missing parameter for argument %s", argument.Data());
260     iResult=-EINVAL;
261   }
262   if (iResult>=0) {
263     iResult=CheckOCDBEntries();
264   }
265   if (iResult>=0) {
266     iResult=DoInit(iNofChildArgs, pArguments);
267   }
268   if (iResult>=0) {
269     fEventCount=0;
270
271     // find out if the component wants to get the steering events
272     // explicitly
273     AliHLTComponentDataTypeList inputDt;
274     GetInputDataTypes(inputDt);
275     bool bRequireSteeringBlocks=false;
276     for (AliHLTComponentDataTypeList::iterator dt=inputDt.begin();
277          dt!=inputDt.end() && !bRequireSteeringBlocks;
278          dt++) {
279       bRequireSteeringBlocks|=MatchExactly(*dt,kAliHLTDataTypeSOR);
280       bRequireSteeringBlocks|=MatchExactly(*dt,kAliHLTDataTypeRunType);
281       bRequireSteeringBlocks|=MatchExactly(*dt,kAliHLTDataTypeEOR);
282       bRequireSteeringBlocks|=MatchExactly(*dt,kAliHLTDataTypeDDL);
283       bRequireSteeringBlocks|=MatchExactly(*dt,kAliHLTDataTypeComponentStatistics);
284     }
285     if (bRequireSteeringBlocks) fFlags|=kRequireSteeringBlocks;
286   }
287   if (pArguments) delete [] pArguments;
288
289 #if defined(HLT_COMPONENT_STATISTICS)
290   // benchmarking stopwatch for the component statistics
291   fpBenchmark=new TStopwatch;
292 #endif // HLT_COMPONENT_STATISTICS
293
294   return iResult;
295 }
296
297 int AliHLTComponent::Deinit()
298 {
299   // see header file for function documentation
300   HLTLogKeyword(fChainId.c_str());
301   int iResult=0;
302   iResult=DoDeinit();
303   if (fpRunDesc) {
304     // TODO: the warning should be kept, but the condition is wrong since the
305     // AliHLTRunDesc is set before the SOR event in the SetRunDescription
306     // method. A couple of state flags should be defined but that is a bit more
307     // work to do. For the moment disable the warning (2009-07-01)
308     // 2009-09-08: now, the info is not cleared in the ProcessEvent, because it
309     // might be needed by components during the event processing.
310     //HLTWarning("did not receive EOR for run %d", fpRunDesc->fRunNo);
311     AliHLTRunDesc* pRunDesc=fpRunDesc;
312     fpRunDesc=NULL;
313     delete pRunDesc;
314   }
315   if (fpCTPData) {
316     delete fpCTPData;
317   }
318   fpCTPData=NULL;
319
320   fEventCount=0;
321   fFlags=0;
322   return iResult;
323 }
324
325 int AliHLTComponent::InitCDB(const char* cdbPath, AliHLTComponentHandler* pHandler)
326 {
327   // see header file for function documentation
328   int iResult=0;
329   HLTInfo("Using CDB: %s", cdbPath);
330   if (pHandler) {
331   // I have to think about separating the library handling from the
332   // component handler. Requiring the component handler here is not
333   // the cleanest solution.
334   // We presume the library already to be loaded, which is the case
335   // because it is loaded in the initialization of the logging functionality
336   //
337   // find the symbol
338   AliHLTMiscInitCDB_t pFunc=(AliHLTMiscInitCDB_t)pHandler->FindSymbol(ALIHLTMISC_LIBRARY, ALIHLTMISC_INIT_CDB);
339   if (pFunc) {
340     TString path;
341     if (cdbPath && cdbPath[0]!=0) {
342       path=cdbPath;
343       // very temporary fix, have to check for other formats
344       if (!path.BeginsWith("local://")) {
345         path="local://";
346         path+=cdbPath;
347       }
348     }
349     if ((iResult=(*pFunc)(path.Data()))>=0) {
350       if (!(fCDBSetRunNoFunc=pHandler->FindSymbol(ALIHLTMISC_LIBRARY, ALIHLTMISC_SET_CDB_RUNNO))) {
351         Message(NULL, kHLTLogWarning, "AliHLTComponent::InitCDB", "init CDB",
352                 "can not find function to set CDB run no");
353       }
354     }
355   } else {
356     Message(NULL, kHLTLogError, "AliHLTComponent::InitCDB", "init CDB",
357             "can not find initialization function");
358     iResult=-ENOSYS;
359   }
360   } else {
361     iResult=-EINVAL;
362   }
363   return iResult;
364 }
365
366 int AliHLTComponent::SetCDBRunNo(int runNo)
367 {
368   // see header file for function documentation
369   if (!fCDBSetRunNoFunc) return 0;
370   return (*((AliHLTMiscSetCDBRunNo_t)fCDBSetRunNoFunc))(runNo);
371 }
372
373 int AliHLTComponent::SetRunDescription(const AliHLTRunDesc* desc, const char* /*runType*/)
374 {
375   // see header file for function documentation
376   if (!desc) return -EINVAL;
377   if (desc->fStructSize!=sizeof(AliHLTRunDesc)) {
378     HLTError("invalid size of RunDesc struct (%ul)", desc->fStructSize);
379     return -EINVAL;
380   }
381
382   if (!fpRunDesc) {
383     fpRunDesc=new AliHLTRunDesc;
384     if (!fpRunDesc) return -ENOMEM;
385     *fpRunDesc=kAliHLTVoidRunDesc;
386   }
387
388   if (fpRunDesc->fRunNo!=kAliHLTVoidRunNo && fpRunDesc->fRunNo!=desc->fRunNo) {
389     HLTWarning("Run description has already been set");
390   }
391   *fpRunDesc=*desc;
392   SetCDBRunNo(fpRunDesc->fRunNo);
393   // TODO: we have to decide about the runType
394   return 0;
395 }
396
397 int AliHLTComponent::SetComponentDescription(const char* desc)
398 {
399   // see header file for function documentation
400   int iResult=0;
401   if (!desc) return 0;
402
403   TString descriptor=desc;
404   TObjArray* pTokens=descriptor.Tokenize(" ");
405   if (pTokens) {
406     for (int i=0; i<pTokens->GetEntriesFast() && iResult>=0; i++) {
407       TString argument=((TObjString*)pTokens->At(i++))->GetString();
408       if (!argument || argument.IsNull()) continue;
409
410       // chainid
411       if (argument.BeginsWith("chainid")) {
412         argument.ReplaceAll("chainid", "");
413         if (argument.BeginsWith("=")) {
414           fChainId=argument.Replace(0,1,"");
415           fChainIdCrc=CalculateChecksum((const AliHLTUInt8_t*)fChainId.c_str(), fChainId.length());
416           HLTDebug("setting component description: chain id %s crc 0x%8x", fChainId.c_str(), fChainIdCrc);
417         } else {
418           fChainId="";
419         }
420       } else {
421         HLTWarning("unknown component description %s", argument.Data());
422       }
423     }
424     delete pTokens;
425   }
426   
427   return iResult;
428 }
429
430 int AliHLTComponent::ConfigureFromArgumentString(int argc, const char** argv)
431 {
432   // see header file for function documentation
433   int iResult=0;
434   vector<const char*> array;
435   TObjArray choppedArguments;
436   TString argument="";
437   int i=0;
438   for (i=0; i<argc && iResult>=0; i++) {
439     argument=argv[i];
440     if (argument.IsWhitespace()) continue;
441
442     // special handling for single component arguments ending with
443     // a sequence of blanks
444     argument.Remove(0, argument.First(' '));
445     if (argument.IsWhitespace()) {
446       array.push_back(argv[i]);
447       continue;
448     }
449
450     // extra blank to insert blank token before leading quotes
451     argument=" ";
452     argument+=argv[i];
453     // insert blank in consecutive quotes to correctly tokenize
454     argument.ReplaceAll("''", "' '");
455     // replace newlines by blanks
456     argument.ReplaceAll("\n", " ");
457     if (argument.IsNull()) continue;
458     TObjArray* pTokensQuote=argument.Tokenize("'");
459     if (pTokensQuote) {
460       if (pTokensQuote->GetEntriesFast()>0) {
461         for (int k=0; k<pTokensQuote->GetEntriesFast(); k++) {
462           argument=((TObjString*)pTokensQuote->At(k))->GetString();
463           if (argument.IsWhitespace()) continue;
464           if (k%2) {
465             // every second entry is enclosed by quotes and thus
466             // one single argument
467             array.push_back(argument.Data());
468           } else {
469     TObjArray* pTokens=argument.Tokenize(" ");
470     if (pTokens) {
471       if (pTokens->GetEntriesFast()>0) {
472         for (int n=0; n<pTokens->GetEntriesFast(); n++) {
473           choppedArguments.AddLast(pTokens->At(n));
474           TString data=((TObjString*)pTokens->At(n))->GetString();
475           if (!data.IsNull()) {
476             array.push_back(data.Data());
477           }
478         }
479         pTokens->SetOwner(kFALSE);
480       }
481       delete pTokens;
482     }
483           }
484         }
485         pTokensQuote->SetOwner(kFALSE);
486       }
487       delete pTokensQuote;
488     }
489   }
490
491   for (i=0; (unsigned)i<array.size() && iResult>=0;) {
492     int result=ScanConfigurationArgument(array.size()-i, &array[i]);
493     if (result==0) {
494       HLTWarning("unknown component argument %s", array[i]);
495       i++;
496     } else if (result>0) {
497       i+=result;
498     } else {
499       iResult=result;
500       if (iResult==-EINVAL) {
501         HLTError("unknown argument %s", array[i]);
502       } else if (iResult==-EPROTO) {
503         HLTError("missing/wrong parameter for argument %s (%s)", array[i], (array.size()>(unsigned)i+1)?array[i+1]:"missing");
504       } else {
505         HLTError("scan of argument %s failed (%d)", array[i], iResult);
506       }
507     }
508   }
509
510   return iResult;
511 }
512
513 int AliHLTComponent::ConfigureFromCDBTObjString(const char* entries)
514 {
515   // see header file for function documentation
516   int iResult=0;
517   TString arguments;
518   TString confEntries=entries;
519   TObjArray* pTokens=confEntries.Tokenize(" ");
520   if (pTokens) {
521     for (int n=0; n<pTokens->GetEntriesFast(); n++) {
522       const char* path=((TObjString*)pTokens->At(n))->GetString().Data();
523       const char* chainId=GetChainId();
524       HLTInfo("configure from entry %s, chain id %s", path, (chainId!=NULL && chainId[0]!=0)?chainId:"<none>");
525       TObject* pOCDBObject = LoadAndExtractOCDBObject(path);
526       if (pOCDBObject) {
527         TObjString* pString=dynamic_cast<TObjString*>(pOCDBObject);
528         if (pString) {
529           HLTInfo("received configuration object string: \'%s\'", pString->GetString().Data());
530           arguments+=pString->GetString().Data();
531           arguments+=" ";
532         } else {
533           HLTError("configuration object \"%s\" has wrong type, required TObjString", path);
534           iResult=-EINVAL;
535         }
536       } else {
537         HLTError("can not fetch object \"%s\" from OCDB", path);
538         iResult=-ENOENT;
539       }
540     }
541     delete pTokens;
542   }
543   if (iResult>=0 && !arguments.IsNull())  {
544     const char* array=arguments.Data();
545     iResult=ConfigureFromArgumentString(1, &array);
546   }
547   return iResult;
548 }
549
550 TObject* AliHLTComponent::LoadAndExtractOCDBObject(const char* path, int version, int subVersion)
551 {
552   // see header file for function documentation
553   AliCDBEntry* pEntry=AliHLTMisc::Instance().LoadOCDBEntry(path, GetRunNo(), version, subVersion);
554   if (!pEntry) return NULL;
555   return AliHLTMisc::Instance().ExtractObject(pEntry);
556 }
557
558 int AliHLTComponent::DoInit( int /*argc*/, const char** /*argv*/)
559 {
560   // default implementation, childs can overload
561   HLTLogKeyword("dummy");
562   return 0;
563 }
564
565 int AliHLTComponent::DoDeinit()
566 {
567   // default implementation, childs can overload
568   HLTLogKeyword("dummy");
569   return 0;
570 }
571
572 int AliHLTComponent::Reconfigure(const char* /*cdbEntry*/, const char* /*chainId*/)
573 {
574   // default implementation, childs can overload
575   HLTLogKeyword("dummy");
576   return 0;
577 }
578
579 int AliHLTComponent::ReadPreprocessorValues(const char* /*modules*/)
580 {
581   // default implementation, childs can overload
582   HLTLogKeyword("dummy");
583   return 0;
584 }
585
586 int AliHLTComponent::ScanConfigurationArgument(int /*argc*/, const char** /*argv*/)
587 {
588   // default implementation, childs can overload
589   HLTLogKeyword("dummy");
590   HLTWarning("The function needs to be implemented by the component");
591   return 0;
592 }
593
594 int AliHLTComponent::StartOfRun()
595 {
596   // default implementation, childs can overload
597   HLTLogKeyword("dummy");
598   return 0;
599 }
600
601 int AliHLTComponent::EndOfRun()
602 {
603   // default implementation, childs can overload
604   HLTLogKeyword("dummy");
605   return 0;
606 }
607
608
609 int AliHLTComponent::GetOutputDataTypes(AliHLTComponentDataTypeList& /*tgtList*/)
610 {
611   // default implementation, childs can overload
612   HLTLogKeyword("dummy");
613   return 0;
614 }
615
616 void AliHLTComponent::GetOCDBObjectDescription( TMap* const /*targetArray*/)
617 {
618   // default implementation, childs can overload
619   HLTLogKeyword("dummy");
620 }
621
622 int AliHLTComponent::CheckOCDBEntries(const TMap* const externList)
623 {
624   // check the availability of the OCDB entry descriptions in the TMap
625   //  key : complete OCDB path of the entry
626   //  value : auxiliary object - short description
627   // if the external map was not provided the function invokes
628   // interface function GetOCDBObjectDescription() to retrieve the list.
629   int iResult=0;
630   if (externList) {
631     iResult=AliHLTMisc::Instance().CheckOCDBEntries(externList);
632   } else {
633     TMap* pMap=new TMap;
634     if (pMap) {
635       pMap->SetOwnerKeyValue(kTRUE);
636       GetOCDBObjectDescription(pMap);
637       iResult=AliHLTMisc::Instance().CheckOCDBEntries(pMap);
638       delete pMap;
639       pMap=NULL;
640     }
641   }
642
643   return iResult;
644 }
645
646 void AliHLTComponent::DataType2Text( const AliHLTComponentDataType& type, char output[kAliHLTComponentDataTypefIDsize+kAliHLTComponentDataTypefOriginSize+2] ) const
647 {
648   // see header file for function documentation
649   memset( output, 0, kAliHLTComponentDataTypefIDsize+kAliHLTComponentDataTypefOriginSize+2 );
650   strncat( output, type.fOrigin, kAliHLTComponentDataTypefOriginSize );
651   strcat( output, ":" );
652   strncat( output, type.fID, kAliHLTComponentDataTypefIDsize );
653 }
654
655 string AliHLTComponent::DataType2Text( const AliHLTComponentDataType& type, int mode)
656 {
657   // see header file for function documentation
658   string out("");
659
660   // 'typeid' 'origin'
661   // aligned to 8 and 4 chars respectively, blocks enclosed in quotes and
662   // separated by blank e.g.
663   // 'DDL_RAW ' 'TPC '
664   if (mode==3) {
665     int i=0;
666     char tmp[8];
667     out+="'";
668     for (i=0; i<kAliHLTComponentDataTypefIDsize; i++) {
669       unsigned char* puc=(unsigned char*)type.fID;
670       if (puc[i]<32)
671         sprintf(tmp, "\\%x", type.fID[i]);
672       else
673         sprintf(tmp, "%c", type.fID[i]);
674       out+=tmp;
675     }
676     out+="' '";
677     for (i=0; i<kAliHLTComponentDataTypefOriginSize; i++) {
678       unsigned char* puc=(unsigned char*)type.fOrigin;
679       if ((puc[i])<32)
680         sprintf(tmp, "\\%x", type.fOrigin[i]);
681       else
682         sprintf(tmp, "%c", type.fOrigin[i]);
683       out+=tmp;
684     }
685     out+="'";
686     return out;
687   }
688
689   // origin typeid as numbers separated by colon e.g.
690   // aligned to 8 and 4 chars respectively, all characters separated by
691   // quotes, e.g.
692   // '84'80'67'32':'68'68'76'95'82'65'87'32'
693   if (mode==2) {
694     int i=0;
695     char tmp[8];
696     for (i=0; i<kAliHLTComponentDataTypefOriginSize; i++) {
697       sprintf(tmp, "'%d", type.fOrigin[i]);
698       out+=tmp;
699     }
700     out+="':'";
701     for (i=0; i<kAliHLTComponentDataTypefIDsize; i++) {
702       sprintf(tmp, "%d'", type.fID[i]);
703       out+=tmp;
704     }
705     return out;
706   }
707
708   // origin typeid separated by colon e.g.
709   // aligned to 8 and 4 chars respectively, all characters separated by
710   // quotes, e.g.
711   // 'T'P'C' ':'D'D'L'_'R'A'W' '
712   if (mode==1) {
713     int i=0;
714     char tmp[8];
715     for (i=0; i<kAliHLTComponentDataTypefOriginSize; i++) {
716       unsigned char* puc=(unsigned char*)type.fOrigin;
717       if ((puc[i])<32)
718         sprintf(tmp, "'\\%x", type.fOrigin[i]);
719       else
720         sprintf(tmp, "'%c", type.fOrigin[i]);
721       out+=tmp;
722     }
723     out+="':'";
724     for (i=0; i<kAliHLTComponentDataTypefIDsize; i++) {
725       unsigned char* puc=(unsigned char*)type.fID;
726       if (puc[i]<32)
727         sprintf(tmp, "\\%x'", type.fID[i]);
728       else
729         sprintf(tmp, "%c'", type.fID[i]);
730       out+=tmp;
731     }
732     return out;
733   }
734
735   // origin typeid
736   // aligned to 8 and 4 chars respectively, separated by colon e.g.
737   // TPC :DDL_RAW 
738   if (type==kAliHLTVoidDataType) {
739     out="VOID:VOID";
740   } else {
741     // some gymnastics in order to avoid a '0' which is part of either or both
742     // ID and origin terminating the whole string. Unfortunately, string doesn't
743     // stop appending at the '0' if the number of elements to append was 
744     // explicitely specified
745     string tmp("");
746     tmp.append(type.fOrigin, kAliHLTComponentDataTypefOriginSize);
747     out.append(tmp.c_str());
748     out.append(":");
749     tmp="";
750     tmp.append(type.fID, kAliHLTComponentDataTypefIDsize);
751     out.append(tmp.c_str());
752   }
753   return out;
754 }
755
756
757 void* AliHLTComponent::AllocMemory( unsigned long size ) 
758 {
759   // see header file for function documentation
760   if (fEnvironment.fAllocMemoryFunc)
761     return (*fEnvironment.fAllocMemoryFunc)(fEnvironment.fParam, size );
762   HLTFatal("no memory allocation handler registered");
763   return NULL;
764 }
765
766 int AliHLTComponent::MakeOutputDataBlockList( const AliHLTComponentBlockDataList& blocks, AliHLTUInt32_t* blockCount,
767                                               AliHLTComponentBlockData** outputBlocks ) 
768 {
769   // see header file for function documentation
770     if ( blockCount==NULL || outputBlocks==NULL )
771         return -EFAULT;
772     AliHLTUInt32_t count = blocks.size();
773     if ( !count )
774         {
775         *blockCount = 0;
776         *outputBlocks = NULL;
777         return 0;
778         }
779     *outputBlocks = reinterpret_cast<AliHLTComponentBlockData*>( AllocMemory( sizeof(AliHLTComponentBlockData)*count ) );
780     if ( !*outputBlocks )
781         return -ENOMEM;
782     for ( unsigned long i = 0; i < count; i++ ) {
783         (*outputBlocks)[i] = blocks[i];
784         if (MatchExactly(blocks[i].fDataType, kAliHLTAnyDataType)) {
785           (*outputBlocks)[i].fDataType=GetOutputDataType();
786           /* data type was set to the output data type by the PubSub AliRoot
787              Wrapper component, if data type of the block was ********:****.
788              Now handled by the component base class in order to have same
789              behavior when running embedded in AliRoot
790           memset((*outputBlocks)[i].fDataType.fID, '*', kAliHLTComponentDataTypefIDsize);
791           memset((*outputBlocks)[i].fDataType.fOrigin, '*', kAliHLTComponentDataTypefOriginSize);
792           */
793         }
794     }
795     *blockCount = count;
796     return 0;
797
798 }
799
800 int AliHLTComponent::GetEventDoneData( unsigned long size, AliHLTComponentEventDoneData** edd ) const
801 {
802   // see header file for function documentation
803   if (fEnvironment.fGetEventDoneDataFunc)
804     return (*fEnvironment.fGetEventDoneDataFunc)(fEnvironment.fParam, fCurrentEvent, size, edd );
805   return -ENOSYS;
806 }
807
808 int AliHLTComponent::ReserveEventDoneData( unsigned long size )
809 {
810   // see header file for function documentation
811   int iResult=0;
812
813   unsigned long capacity=fEventDoneDataSize;
814   if (fEventDoneData) capacity-=sizeof(AliHLTComponentEventDoneData)+fEventDoneData->fDataSize;
815   if (size>capacity) {
816     unsigned long newSize=sizeof(AliHLTComponentEventDoneData)+size+(fEventDoneDataSize-capacity);
817     AliHLTComponentEventDoneData* newEDD = reinterpret_cast<AliHLTComponentEventDoneData*>( new AliHLTUInt8_t[newSize] );
818     if (!newEDD)
819       return -ENOMEM;
820     newEDD->fStructSize = sizeof(AliHLTComponentEventDoneData);
821     newEDD->fDataSize = 0;
822     newEDD->fData = reinterpret_cast<AliHLTUInt8_t*>(newEDD)+newEDD->fStructSize;
823     if (fEventDoneData) {
824       memcpy( newEDD->fData, fEventDoneData->fData, fEventDoneData->fDataSize );
825       newEDD->fDataSize = fEventDoneData->fDataSize;
826       delete [] reinterpret_cast<AliHLTUInt8_t*>( fEventDoneData );
827     }
828     fEventDoneData = newEDD;
829     fEventDoneDataSize = newSize;
830   }
831   return iResult;
832
833 }
834
835 int AliHLTComponent::PushEventDoneData( AliHLTUInt32_t eddDataWord )
836 {
837   // see header file for function documentation
838   if (!fEventDoneData)
839     return -ENOMEM;
840   if (fEventDoneData->fDataSize+sizeof(AliHLTUInt32_t)>fEventDoneDataSize)
841     return -ENOSPC;
842   *reinterpret_cast<AliHLTUInt32_t*>((reinterpret_cast<AliHLTUInt8_t*>(fEventDoneData->fData)+fEventDoneData->fDataSize)) = eddDataWord;
843   fEventDoneData->fDataSize += sizeof(AliHLTUInt32_t);
844   return 0;
845 }
846
847 void AliHLTComponent::ReleaseEventDoneData()
848 {
849    // see header file for function documentation
850  if (fEventDoneData)
851     delete [] reinterpret_cast<AliHLTUInt8_t*>( fEventDoneData );
852   fEventDoneData = NULL;
853   fEventDoneDataSize = 0;
854 }
855
856
857 int AliHLTComponent::FindMatchingDataTypes(AliHLTComponent* pConsumer, AliHLTComponentDataTypeList* tgtList) 
858 {
859   // see header file for function documentation
860   int iResult=0;
861   if (pConsumer) {
862     AliHLTComponentDataTypeList itypes;
863     AliHLTComponentDataTypeList otypes;
864     otypes.push_back(GetOutputDataType());
865     if (MatchExactly(otypes[0],kAliHLTMultipleDataType)) {
866       otypes.clear();
867       int count=0;
868       if ((count=GetOutputDataTypes(otypes))>0) {
869       } else if (GetComponentType()!=kSink) {
870         HLTWarning("component %s indicates multiple output data types but GetOutputDataTypes returns %d", GetComponentID(), count);
871       }
872     }
873     pConsumer->GetInputDataTypes(itypes);
874     AliHLTComponentDataTypeList::iterator otype=otypes.begin();
875     for (;otype!=otypes.end();otype++) {
876       //PrintDataTypeContent((*otype), "publisher \'%s\'");
877       if ((*otype)==(kAliHLTAnyDataType|kAliHLTDataOriginPrivate)) {
878         if (tgtList) tgtList->push_back(*otype);
879         iResult++;
880         continue;
881       }
882       
883       AliHLTComponentDataTypeList::iterator itype=itypes.begin();
884       for ( ; itype!=itypes.end() && (*itype)!=(*otype) ; itype++) {/* empty body */};
885       //if (itype!=itypes.end()) PrintDataTypeContent(*itype, "consumer \'%s\'");
886       if (itype!=itypes.end()) {
887         if (tgtList) tgtList->push_back(*otype);
888         iResult++;
889       }
890     }
891   } else {
892     iResult=-EINVAL;
893   }
894   return iResult;
895 }
896
897 void AliHLTComponent::PrintDataTypeContent(AliHLTComponentDataType& dt, const char* format)
898 {
899   // see header file for function documentation
900   const char* fmt="\'%s\'";
901   if (format) fmt=format;
902   AliHLTLogging::Message(NULL, kHLTLogNone, NULL , NULL, Form(fmt, (DataType2Text(dt)).c_str()));
903   AliHLTLogging::Message(NULL, kHLTLogNone, NULL , NULL, 
904                          Form("%x %x %x %x %x %x %x %x : %x %x %x %x", 
905                               dt.fID[0],
906                               dt.fID[1],
907                               dt.fID[2],
908                               dt.fID[3],
909                               dt.fID[4],
910                               dt.fID[5],
911                               dt.fID[6],
912                               dt.fID[7],
913                               dt.fOrigin[0],
914                               dt.fOrigin[1],
915                               dt.fOrigin[2],
916                               dt.fOrigin[3]));
917 }
918
919 void AliHLTComponent::FillBlockData( AliHLTComponentBlockData& blockData )
920 {
921   // see header file for function documentation
922   blockData.fStructSize = sizeof(blockData);
923   FillShmData( blockData.fShmKey );
924   blockData.fOffset = ~(AliHLTUInt32_t)0;
925   blockData.fPtr = NULL;
926   blockData.fSize = 0;
927   FillDataType( blockData.fDataType );
928   blockData.fSpecification = kAliHLTVoidDataSpec;
929 }
930
931 void AliHLTComponent::FillShmData( AliHLTComponentShmData& shmData )
932 {
933   // see header file for function documentation
934   shmData.fStructSize = sizeof(shmData);
935   shmData.fShmType = gkAliHLTComponentInvalidShmType;
936   shmData.fShmID = gkAliHLTComponentInvalidShmID;
937 }
938
939 void AliHLTComponent::FillDataType( AliHLTComponentDataType& dataType )
940 {
941   // see header file for function documentation
942   dataType=kAliHLTAnyDataType;
943 }
944
945 void AliHLTComponent::CopyDataType(AliHLTComponentDataType& tgtdt, const AliHLTComponentDataType& srcdt) 
946 {
947   // see header file for function documentation
948   memcpy(&tgtdt.fID[0], &srcdt.fID[0], kAliHLTComponentDataTypefIDsize);
949   memcpy(&tgtdt.fOrigin[0], &srcdt.fOrigin[0], kAliHLTComponentDataTypefOriginSize);
950 }
951
952 void AliHLTComponent::SetDataType(AliHLTComponentDataType& tgtdt, const char* id, const char* origin) 
953 {
954   // see header file for function documentation
955   tgtdt.fStructSize=sizeof(AliHLTComponentDataType);
956   if (id) {
957     memset(&tgtdt.fID[0], 0, kAliHLTComponentDataTypefIDsize);
958     strncpy(&tgtdt.fID[0], id, strlen(id)<(size_t)kAliHLTComponentDataTypefIDsize?strlen(id):kAliHLTComponentDataTypefIDsize);
959   }
960   if (origin) {
961     memset(&tgtdt.fOrigin[0], 0, kAliHLTComponentDataTypefOriginSize);
962     strncpy(&tgtdt.fOrigin[0], origin, strlen(origin)<(size_t)kAliHLTComponentDataTypefOriginSize?strlen(origin):kAliHLTComponentDataTypefOriginSize);
963   }
964 }
965
966 void AliHLTComponent::SetDataType(AliHLTComponentDataType& dt, AliHLTUInt64_t id, AliHLTUInt32_t origin)
967 {
968   // see header file for function documentation
969   dt.fStructSize=sizeof(AliHLTComponentDataType);
970   assert(kAliHLTComponentDataTypefIDsize==sizeof(id));
971   assert(kAliHLTComponentDataTypefOriginSize==sizeof(origin));
972   memcpy(&dt.fID, &id, kAliHLTComponentDataTypefIDsize);
973   memcpy(&dt.fOrigin, &origin, kAliHLTComponentDataTypefOriginSize);
974 }
975
976 void AliHLTComponent::FillEventData(AliHLTComponentEventData& evtData)
977 {
978   // see header file for function documentation
979   memset(&evtData, 0, sizeof(AliHLTComponentEventData));
980   evtData.fStructSize=sizeof(AliHLTComponentEventData);
981   evtData.fEventID=kAliHLTVoidEventID;
982 }
983
984 void AliHLTComponent::PrintComponentDataTypeInfo(const AliHLTComponentDataType& dt) 
985 {
986   // see header file for function documentation
987   TString msg;
988   msg.Form("AliHLTComponentDataType(%d): ID=\"", dt.fStructSize);
989   for ( int i = 0; i < kAliHLTComponentDataTypefIDsize; i++ ) {
990    if (dt.fID[i]!=0) msg+=dt.fID[i];
991    else msg+="\\0";
992   }
993   msg+="\" Origin=\"";
994   for ( int i = 0; i < kAliHLTComponentDataTypefOriginSize; i++ ) {
995    if (dt.fOrigin[i]!=0) msg+=dt.fOrigin[i];
996    else msg+="\\0";
997   }
998   msg+="\"";
999   AliHLTLogging::Message(NULL, kHLTLogNone, NULL , NULL, msg.Data());
1000 }
1001
1002 int AliHLTComponent::GetEventCount() const
1003 {
1004   // see header file for function documentation
1005   return fEventCount;
1006 }
1007
1008 int AliHLTComponent::IncrementEventCounter()
1009 {
1010   // see header file for function documentation
1011   if (fEventCount>=0) fEventCount++;
1012   return fEventCount;
1013 }
1014
1015 int AliHLTComponent::GetNumberOfInputBlocks() const
1016 {
1017   // see header file for function documentation
1018   if (fpInputBlocks!=NULL) {
1019     return fCurrentEventData.fBlockCnt;
1020   }
1021   return 0;
1022 }
1023
1024 AliHLTEventID_t AliHLTComponent::GetEventId() const
1025 {
1026   // see header file for function documentation
1027   if (fpInputBlocks!=NULL) {
1028     return fCurrentEventData.fEventID;
1029   }
1030   return 0;
1031 }
1032
1033 const TObject* AliHLTComponent::GetFirstInputObject(const AliHLTComponentDataType& dt,
1034                                                     const char* classname,
1035                                                     int bForce)
1036 {
1037   // see header file for function documentation
1038   ALIHLTCOMPONENT_BASE_STOPWATCH();
1039   fSearchDataType=dt;
1040   if (classname) fClassName=classname;
1041   else fClassName.clear();
1042   int idx=FindInputBlock(fSearchDataType, 0, 1);
1043   TObject* pObj=NULL;
1044   if (idx>=0) {
1045     HLTDebug("found block %d when searching for data type %s", idx, DataType2Text(dt).c_str());
1046     if ((pObj=GetInputObject(idx, fClassName.c_str(), bForce))!=NULL) {
1047       fCurrentInputBlock=idx;
1048     } else {
1049     }
1050   }
1051   return pObj;
1052 }
1053
1054 const TObject* AliHLTComponent::GetFirstInputObject(const char* dtID, 
1055                                                     const char* dtOrigin,
1056                                                     const char* classname,
1057                                                     int         bForce)
1058 {
1059   // see header file for function documentation
1060   ALIHLTCOMPONENT_BASE_STOPWATCH();
1061   AliHLTComponentDataType dt;
1062   SetDataType(dt, dtID, dtOrigin);
1063   return GetFirstInputObject(dt, classname, bForce);
1064 }
1065
1066 const TObject* AliHLTComponent::GetNextInputObject(int bForce)
1067 {
1068   // see header file for function documentation
1069   ALIHLTCOMPONENT_BASE_STOPWATCH();
1070   int idx=FindInputBlock(fSearchDataType, fCurrentInputBlock+1, 1);
1071   //HLTDebug("found block %d when searching for data type %s", idx, DataType2Text(fSearchDataType).c_str());
1072   TObject* pObj=NULL;
1073   if (idx>=0) {
1074     if ((pObj=GetInputObject(idx, fClassName.c_str(), bForce))!=NULL) {
1075       fCurrentInputBlock=idx;
1076     }
1077   }
1078   return pObj;
1079 }
1080
1081 int AliHLTComponent::FindInputBlock(const AliHLTComponentDataType& dt, int startIdx, int bObject) const
1082 {
1083   // see header file for function documentation
1084   int iResult=-ENOENT;
1085   if (fpInputBlocks!=NULL) {
1086     int idx=startIdx<0?0:startIdx;
1087     for ( ; (UInt_t)idx<fCurrentEventData.fBlockCnt && iResult==-ENOENT; idx++) {
1088       if (dt!=fpInputBlocks[idx].fDataType) continue;
1089
1090       if (bObject!=0) {
1091         if (fpInputBlocks[idx].fPtr==NULL) continue;
1092         AliHLTUInt32_t firstWord=*((AliHLTUInt32_t*)fpInputBlocks[idx].fPtr);
1093         if (firstWord!=fpInputBlocks[idx].fSize-sizeof(AliHLTUInt32_t)) continue;
1094       }
1095       iResult=idx;
1096     }
1097   }
1098   return iResult;
1099 }
1100
1101 TObject* AliHLTComponent::CreateInputObject(int idx, int bForce)
1102 {
1103   // see header file for function documentation
1104   TObject* pObj=NULL;
1105   if (fpInputBlocks!=NULL) {
1106     if ((UInt_t)idx<fCurrentEventData.fBlockCnt) {
1107       if (fpInputBlocks[idx].fPtr) {
1108         AliHLTUInt32_t firstWord=*((AliHLTUInt32_t*)fpInputBlocks[idx].fPtr);
1109         if (firstWord==fpInputBlocks[idx].fSize-sizeof(AliHLTUInt32_t)) {
1110           HLTDebug("create object from block %d size %d", idx, fpInputBlocks[idx].fSize);
1111           AliHLTMessage msg(fpInputBlocks[idx].fPtr, fpInputBlocks[idx].fSize);
1112           TClass* objclass=msg.GetClass();
1113           pObj=msg.ReadObject(objclass);
1114           if (pObj && objclass) {
1115             HLTDebug("object %p type %s created", pObj, objclass->GetName());
1116           } else {
1117           }
1118           //} else {
1119         } else if (bForce!=0) {
1120           HLTError("size mismatch: block size %d, indicated %d", fpInputBlocks[idx].fSize, firstWord+sizeof(AliHLTUInt32_t));
1121         }
1122       } else {
1123         HLTFatal("block descriptor empty");
1124       }
1125     } else {
1126       HLTError("index %d out of range %d", idx, fCurrentEventData.fBlockCnt);
1127     }
1128   } else {
1129     HLTError("no input blocks available");
1130   }
1131   
1132   return pObj;
1133 }
1134
1135 TObject* AliHLTComponent::GetInputObject(int idx, const char* /*classname*/, int bForce)
1136 {
1137   // see header file for function documentation
1138   if (fpInputObjects==NULL) {
1139     fpInputObjects=new TObjArray(fCurrentEventData.fBlockCnt);
1140   }
1141   TObject* pObj=NULL;
1142   if (fpInputObjects) {
1143     pObj=fpInputObjects->At(idx);
1144     if (pObj==NULL) {
1145       pObj=CreateInputObject(idx, bForce);
1146       if (pObj) {
1147         fpInputObjects->AddAt(pObj, idx);
1148       }
1149     }
1150   } else {
1151     HLTFatal("memory allocation failed: TObjArray of size %d", fCurrentEventData.fBlockCnt);
1152   }
1153   return pObj;
1154 }
1155
1156 int AliHLTComponent::CleanupInputObjects()
1157 {
1158   // see header file for function documentation
1159   if (!fpInputObjects) return 0;
1160   TObjArray* array=fpInputObjects;
1161   fpInputObjects=NULL;
1162   for (int i=0; i<array->GetEntriesFast(); i++) {
1163     TObject* pObj=array->At(i);
1164     // grrr, garbage collection strikes back: When read via AliHLTMessage
1165     // (CreateInputObject), and written to a TFile afterwards, the
1166     // TFile::Close calls ROOOT's garbage collection. No clue why the
1167     // object ended up in the key list and needs to be deleted
1168     //
1169     // Matthias 09.11.2008 follow up
1170     // This approach doesn't actually work in all cases: the object table
1171     // can be switched off globally, the flag needs to be checked here as
1172     // well in order to avoid memory leaks.
1173     // This means we have to find another solution for the problem if it
1174     // pops up again.
1175     if (pObj &&
1176         (!TObject::GetObjectStat() || gObjectTable->PtrIsValid(pObj))) {
1177       delete pObj;
1178     }
1179   }
1180   delete array;
1181   return 0;
1182 }
1183
1184 AliHLTComponentDataType AliHLTComponent::GetDataType(const TObject* pObject)
1185 {
1186   // see header file for function documentation
1187   ALIHLTCOMPONENT_BASE_STOPWATCH();
1188   AliHLTComponentDataType dt=kAliHLTVoidDataType;
1189   int idx=fCurrentInputBlock;
1190   if (pObject) {
1191     if (fpInputObjects==NULL || (idx=fpInputObjects->IndexOf(pObject))>=0) {
1192     } else {
1193       HLTError("unknown object %p", pObject);
1194     }
1195   }
1196   if (idx>=0) {
1197     if ((UInt_t)idx<fCurrentEventData.fBlockCnt) {
1198       dt=fpInputBlocks[idx].fDataType;
1199     } else {
1200       HLTFatal("severe internal error, index out of range");
1201     }
1202   }
1203   return dt;
1204 }
1205
1206 AliHLTUInt32_t AliHLTComponent::GetSpecification(const TObject* pObject)
1207 {
1208   // see header file for function documentation
1209   ALIHLTCOMPONENT_BASE_STOPWATCH();
1210   AliHLTUInt32_t iSpec=kAliHLTVoidDataSpec;
1211   int idx=fCurrentInputBlock;
1212   if (pObject) {
1213     if (fpInputObjects==NULL || (idx=fpInputObjects->IndexOf(pObject))>=0) {
1214     } else {
1215       HLTError("unknown object %p", pObject);
1216     }
1217   }
1218   if (idx>=0) {
1219     if ((UInt_t)idx<fCurrentEventData.fBlockCnt) {
1220       iSpec=fpInputBlocks[idx].fSpecification;
1221     } else {
1222       HLTFatal("severe internal error, index out of range");
1223     }
1224   }
1225   return iSpec;
1226 }
1227
1228 int AliHLTComponent::Forward(const TObject* pObject)
1229 {
1230   // see header file for function documentation
1231   int iResult=0;
1232   int idx=fCurrentInputBlock;
1233   if (pObject) {
1234     if (fpInputObjects==NULL || (idx=fpInputObjects->IndexOf(pObject))>=0) {
1235     } else {
1236       HLTError("unknown object %p", pObject);
1237       iResult=-ENOENT;
1238     }
1239   }
1240   if (idx>=0) {
1241     fOutputBlocks.push_back(fpInputBlocks[idx]);
1242   }
1243   return iResult;
1244 }
1245
1246 int AliHLTComponent::Forward(const AliHLTComponentBlockData* pBlock)
1247 {
1248   // see header file for function documentation
1249   int iResult=0;
1250   int idx=fCurrentInputBlock;
1251   if (pBlock) {
1252     if ((idx=FindInputBlock(pBlock))>=0) {
1253     } else {
1254       HLTError("unknown Block %p", pBlock);
1255       iResult=-ENOENT;      
1256     }
1257   }
1258   if (idx>=0) {
1259     // check for fpInputBlocks pointer done in FindInputBlock
1260     fOutputBlocks.push_back(fpInputBlocks[idx]);
1261   }
1262   return iResult;
1263 }
1264
1265 const AliHLTComponentBlockData* AliHLTComponent::GetFirstInputBlock(const AliHLTComponentDataType& dt)
1266 {
1267   // see header file for function documentation
1268   ALIHLTCOMPONENT_BASE_STOPWATCH();
1269   fSearchDataType=dt;
1270   fClassName.clear();
1271   int idx=FindInputBlock(fSearchDataType, 0);
1272   const AliHLTComponentBlockData* pBlock=NULL;
1273   if (idx>=0) {
1274     // check for fpInputBlocks pointer done in FindInputBlock
1275     pBlock=&fpInputBlocks[idx];
1276     fCurrentInputBlock=idx;
1277   }
1278   return pBlock;
1279 }
1280
1281 const AliHLTComponentBlockData* AliHLTComponent::GetFirstInputBlock(const char* dtID, 
1282                                                                     const char* dtOrigin)
1283 {
1284   // see header file for function documentation
1285   ALIHLTCOMPONENT_BASE_STOPWATCH();
1286   AliHLTComponentDataType dt;
1287   SetDataType(dt, dtID, dtOrigin);
1288   return GetFirstInputBlock(dt);
1289 }
1290
1291 const AliHLTComponentBlockData* AliHLTComponent::GetInputBlock(int index) const
1292 {
1293   // see header file for function documentation
1294   ALIHLTCOMPONENT_BASE_STOPWATCH();
1295   assert( 0 <= index and index < (int)fCurrentEventData.fBlockCnt );
1296   return &fpInputBlocks[index];
1297 }
1298
1299 const AliHLTComponentBlockData* AliHLTComponent::GetNextInputBlock()
1300 {
1301   // see header file for function documentation
1302   ALIHLTCOMPONENT_BASE_STOPWATCH();
1303   int idx=FindInputBlock(fSearchDataType, fCurrentInputBlock+1);
1304   const AliHLTComponentBlockData* pBlock=NULL;
1305   if (idx>=0) {
1306     // check for fpInputBlocks pointer done in FindInputBlock
1307     pBlock=&fpInputBlocks[idx];
1308     fCurrentInputBlock=idx;
1309   }
1310   return pBlock;
1311 }
1312
1313 int AliHLTComponent::FindInputBlock(const AliHLTComponentBlockData* pBlock) const
1314 {
1315   // see header file for function documentation
1316   int iResult=-ENOENT;
1317   if (fpInputBlocks!=NULL) {
1318     if (pBlock) {
1319       if (pBlock>=fpInputBlocks && pBlock<fpInputBlocks+fCurrentEventData.fBlockCnt) {
1320         iResult=(int)(pBlock-fpInputBlocks);
1321       }
1322     } else {
1323       iResult=-EINVAL;
1324     }
1325   }
1326   return iResult;
1327 }
1328
1329 AliHLTUInt32_t AliHLTComponent::GetSpecification(const AliHLTComponentBlockData* pBlock)
1330 {
1331   // see header file for function documentation
1332   ALIHLTCOMPONENT_BASE_STOPWATCH();
1333   AliHLTUInt32_t iSpec=kAliHLTVoidDataSpec;
1334   int idx=fCurrentInputBlock;
1335   if (pBlock) {
1336     if (fpInputObjects==NULL || (idx=FindInputBlock(pBlock))>=0) {
1337     } else {
1338       HLTError("unknown Block %p", pBlock);
1339     }
1340   }
1341   if (idx>=0) {
1342     // check for fpInputBlocks pointer done in FindInputBlock
1343     iSpec=fpInputBlocks[idx].fSpecification;
1344   }
1345   return iSpec;
1346 }
1347
1348 int AliHLTComponent::PushBack(TObject* pObject, const AliHLTComponentDataType& dt, AliHLTUInt32_t spec, 
1349                               void* pHeader, int headerSize)
1350 {
1351   // see header file for function documentation
1352   ALIHLTCOMPONENT_BASE_STOPWATCH();
1353   int iResult=0;
1354   fLastObjectSize=0;
1355   if (fPushbackPeriod>0) {
1356     // suppress the output
1357     TDatime time;
1358     if (fLastPushBackTime<0 || (int)time.Get()-fLastPushBackTime<fPushbackPeriod) return 0;
1359   }
1360   if (pObject) {
1361     AliHLTMessage msg(kMESS_OBJECT);
1362     msg.SetCompressionLevel(fCompressionLevel);
1363     msg.WriteObject(pObject);
1364     Int_t iMsgLength=msg.Length();
1365     if (iMsgLength>0) {
1366       // Matthias Sep 2008
1367       // NOTE: AliHLTMessage does implement it's own SetLength method
1368       // which is not architecture independent. The original SetLength
1369       // stores the size always in network byte order.
1370       // I'm trying to remember the rational for that, might be that
1371       // it was just some lack of knowledge. Want to change this, but
1372       // has to be done carefullt to be backward compatible.
1373       msg.SetLength(); // sets the length to the first (reserved) word
1374
1375       // does nothing if the level is 0
1376       msg.Compress();
1377
1378       char *mbuf = msg.Buffer();
1379       if (msg.CompBuffer()) {
1380         msg.SetLength(); // set once more to have to byte order
1381         mbuf = msg.CompBuffer();
1382         iMsgLength = msg.CompLength();
1383       }
1384       assert(mbuf!=NULL);
1385       iResult=InsertOutputBlock(mbuf, iMsgLength, dt, spec, pHeader, headerSize);
1386       if (iResult>=0) {
1387         HLTDebug("object %s (%p) size %d compression %d inserted to output", pObject->ClassName(), pObject, iMsgLength, msg.GetCompressionLevel());
1388       }
1389       fLastObjectSize=iMsgLength;
1390     } else {
1391       HLTError("object serialization failed for object %p", pObject);
1392       iResult=-ENOMSG;
1393     }
1394   } else {
1395     iResult=-EINVAL;
1396   }
1397   return iResult;
1398 }
1399
1400 int AliHLTComponent::PushBack(TObject* pObject, const char* dtID, const char* dtOrigin, AliHLTUInt32_t spec,
1401                               void* pHeader, int headerSize)
1402 {
1403   // see header file for function documentation
1404   ALIHLTCOMPONENT_BASE_STOPWATCH();
1405   AliHLTComponentDataType dt;
1406   SetDataType(dt, dtID, dtOrigin);
1407   return PushBack(pObject, dt, spec, pHeader, headerSize);
1408 }
1409
1410 int AliHLTComponent::PushBack(const void* pBuffer, int iSize, const AliHLTComponentDataType& dt, AliHLTUInt32_t spec,
1411                               const void* pHeader, int headerSize)
1412 {
1413   // see header file for function documentation
1414   ALIHLTCOMPONENT_BASE_STOPWATCH();
1415   if (fPushbackPeriod>0) {
1416     // suppress the output
1417     TDatime time;
1418     if (fLastPushBackTime<0 || (int)time.Get()-fLastPushBackTime<fPushbackPeriod) return 0;
1419   }
1420
1421   return InsertOutputBlock(pBuffer, iSize, dt, spec, pHeader, headerSize);
1422 }
1423
1424 int AliHLTComponent::PushBack(const void* pBuffer, int iSize, const char* dtID, const char* dtOrigin, AliHLTUInt32_t spec,
1425                               const void* pHeader, int headerSize)
1426 {
1427   // see header file for function documentation
1428   ALIHLTCOMPONENT_BASE_STOPWATCH();
1429   AliHLTComponentDataType dt;
1430   SetDataType(dt, dtID, dtOrigin);
1431   return PushBack(pBuffer, iSize, dt, spec, pHeader, headerSize);
1432 }
1433
1434 int AliHLTComponent::InsertOutputBlock(const void* pBuffer, int iBufferSize, const AliHLTComponentDataType& dt, AliHLTUInt32_t spec,
1435                                        const void* pHeader, int iHeaderSize)
1436 {
1437   // see header file for function documentation
1438   int iResult=0;
1439   int iBlkSize = iBufferSize + iHeaderSize;
1440
1441   if ((pBuffer!=NULL && iBufferSize>0) || (pHeader!=NULL && iHeaderSize>0)) {
1442     if (fpOutputBuffer && iBlkSize<=(int)(fOutputBufferSize-fOutputBufferFilled)) {
1443       AliHLTUInt8_t* pTgt=fpOutputBuffer+fOutputBufferFilled;
1444
1445       // copy header if provided but skip if the header is the target location
1446       // in that case it has already been copied
1447       if (pHeader!=NULL && pHeader!=pTgt) {
1448         memcpy(pTgt, pHeader, iHeaderSize);
1449       }
1450
1451       pTgt += (AliHLTUInt8_t) iHeaderSize;
1452
1453       // copy buffer if provided but skip if buffer is the target location
1454       // in that case it has already been copied
1455       if (pBuffer!=NULL && pBuffer!=pTgt) {
1456         memcpy(pTgt, pBuffer, iBufferSize);
1457         
1458         //AliHLTUInt32_t firstWord=*((AliHLTUInt32_t*)pBuffer); 
1459         //HLTDebug("copy %d bytes from %p to output buffer %p, first word %#x", iBufferSize, pBuffer, pTgt, firstWord);
1460       }
1461       //HLTDebug("buffer inserted to output: size %d data type %s spec %#x", iBlkSize, DataType2Text(dt).c_str(), spec);
1462     } else {
1463       if (fpOutputBuffer) {
1464         HLTError("too little space in output buffer: %d of %d, required %d", fOutputBufferSize-fOutputBufferFilled, fOutputBufferSize, iBlkSize);
1465       } else {
1466         HLTError("output buffer not available");
1467       }
1468       iResult=-ENOSPC;
1469     }
1470   }
1471   if (iResult>=0) {
1472     AliHLTComponentBlockData bd;
1473     FillBlockData( bd );
1474     bd.fOffset        = fOutputBufferFilled;
1475     bd.fSize          = iBlkSize;
1476     bd.fDataType      = dt;
1477     bd.fSpecification = spec;
1478     fOutputBlocks.push_back( bd );
1479     fOutputBufferFilled+=bd.fSize;
1480   }
1481
1482   return iResult;
1483 }
1484
1485 int AliHLTComponent::EstimateObjectSize(TObject* pObject) const
1486 {
1487   // see header file for function documentation
1488   if (!pObject) return 0;
1489
1490   AliHLTMessage msg(kMESS_OBJECT);
1491   msg.WriteObject(pObject);
1492   return msg.Length();  
1493 }
1494
1495 AliHLTMemoryFile* AliHLTComponent::CreateMemoryFile(int capacity, const char* dtID,
1496                                                     const char* dtOrigin,
1497                                                     AliHLTUInt32_t spec)
1498 {
1499   // see header file for function documentation
1500   ALIHLTCOMPONENT_BASE_STOPWATCH();
1501   AliHLTComponentDataType dt;
1502   SetDataType(dt, dtID, dtOrigin);
1503   return CreateMemoryFile(capacity, dt, spec);
1504 }
1505
1506 AliHLTMemoryFile* AliHLTComponent::CreateMemoryFile(int capacity,
1507                                                     const AliHLTComponentDataType& dt,
1508                                                     AliHLTUInt32_t spec)
1509 {
1510   // see header file for function documentation
1511   ALIHLTCOMPONENT_BASE_STOPWATCH();
1512   AliHLTMemoryFile* pFile=NULL;
1513   if (capacity>=0 && static_cast<unsigned int>(capacity)<=fOutputBufferSize-fOutputBufferFilled){
1514     AliHLTUInt8_t* pTgt=fpOutputBuffer+fOutputBufferFilled;
1515     pFile=new AliHLTMemoryFile((char*)pTgt, capacity);
1516     if (pFile) {
1517       unsigned int nofBlocks=fOutputBlocks.size();
1518       if (nofBlocks+1>fMemFiles.size()) {
1519         fMemFiles.resize(nofBlocks+1, NULL);
1520       }
1521       if (nofBlocks<fMemFiles.size()) {
1522         fMemFiles[nofBlocks]=pFile;
1523         AliHLTComponentBlockData bd;
1524         FillBlockData( bd );
1525         bd.fOffset        = fOutputBufferFilled;
1526         bd.fSize          = capacity;
1527         bd.fDataType      = dt;
1528         bd.fSpecification = spec;
1529         fOutputBufferFilled+=bd.fSize;
1530         fOutputBlocks.push_back( bd );
1531       } else {
1532         HLTError("can not allocate/grow object array");
1533         pFile->CloseMemoryFile(0);
1534         delete pFile;
1535         pFile=NULL;
1536       }
1537     }
1538   } else {
1539     HLTError("can not create memory file of size %d (%d available)", capacity, fOutputBufferSize-fOutputBufferFilled);
1540   }
1541   return pFile;
1542 }
1543
1544 AliHLTMemoryFile* AliHLTComponent::CreateMemoryFile(const char* dtID,
1545                                                     const char* dtOrigin,
1546                                                     AliHLTUInt32_t spec,
1547                                                     float capacity)
1548 {
1549   // see header file for function documentation
1550   ALIHLTCOMPONENT_BASE_STOPWATCH();
1551   AliHLTComponentDataType dt;
1552   SetDataType(dt, dtID, dtOrigin);
1553   int size=fOutputBufferSize-fOutputBufferFilled;
1554   if (capacity<0 || capacity>1.0) {
1555     HLTError("invalid parameter: capacity %f", capacity);
1556     return NULL;
1557   }
1558   size=(int)(size*capacity);
1559   return CreateMemoryFile(size, dt, spec);
1560 }
1561
1562 AliHLTMemoryFile* AliHLTComponent::CreateMemoryFile(const AliHLTComponentDataType& dt,
1563                                                     AliHLTUInt32_t spec,
1564                                                     float capacity)
1565 {
1566   // see header file for function documentation
1567   ALIHLTCOMPONENT_BASE_STOPWATCH();
1568   int size=fOutputBufferSize-fOutputBufferFilled;
1569   if (capacity<0 || capacity>1.0) {
1570     HLTError("invalid parameter: capacity %f", capacity);
1571     return NULL;
1572   }
1573   size=(int)(size*capacity);
1574   return CreateMemoryFile(size, dt, spec);
1575 }
1576
1577 int AliHLTComponent::Write(AliHLTMemoryFile* pFile, const TObject* pObject,
1578                            const char* key, int option)
1579 {
1580   // see header file for function documentation
1581   int iResult=0;
1582   if (pFile && pObject) {
1583     pFile->cd();
1584     iResult=pObject->Write(key, option);
1585     if (iResult>0) {
1586       // success
1587     } else {
1588       iResult=-pFile->GetErrno();
1589       if (iResult==-ENOSPC) {
1590         HLTError("error writing memory file, buffer too small");
1591       }
1592     }
1593   } else {
1594     iResult=-EINVAL;
1595   }
1596   return iResult;
1597 }
1598
1599 int AliHLTComponent::CloseMemoryFile(AliHLTMemoryFile* pFile)
1600 {
1601   // see header file for function documentation
1602   int iResult=0;
1603   if (pFile) {
1604     AliHLTMemoryFilePList::iterator element=fMemFiles.begin();
1605     int i=0;
1606     while (element!=fMemFiles.end() && iResult>=0) {
1607       if (*element && *element==pFile) {
1608         iResult=pFile->CloseMemoryFile();
1609         
1610         // sync memory files and descriptors
1611         if (iResult>=0) {
1612           fOutputBlocks[i].fSize=(*element)->GetSize()+(*element)->GetHeaderSize();
1613         }
1614         delete *element;
1615         *element=NULL;
1616         return iResult;
1617       }
1618       element++; i++;
1619     }
1620     HLTError("can not find memory file %p", pFile);
1621     iResult=-ENOENT;
1622   } else {
1623     iResult=-EINVAL;
1624   }
1625   return iResult;
1626 }
1627
1628 int AliHLTComponent::CreateEventDoneData(AliHLTComponentEventDoneData edd)
1629 {
1630   // see header file for function documentation
1631   int iResult=0;
1632
1633   AliHLTComponentEventDoneData* newEDD = NULL;
1634   
1635   unsigned long newSize=edd.fDataSize;
1636   if (fEventDoneData)
1637     newSize += fEventDoneData->fDataSize;
1638
1639   if (newSize>fEventDoneDataSize) {
1640     newEDD = reinterpret_cast<AliHLTComponentEventDoneData*>( new AliHLTUInt8_t[ sizeof(AliHLTComponentEventDoneData)+newSize ] );
1641     if (!newEDD)
1642       return -ENOMEM;
1643     newEDD->fStructSize = sizeof(AliHLTComponentEventDoneData);
1644     newEDD->fDataSize = newSize;
1645     newEDD->fData = reinterpret_cast<AliHLTUInt8_t*>(newEDD)+newEDD->fStructSize;
1646     unsigned long long offset = 0;
1647     if (fEventDoneData) {
1648       memcpy( newEDD->fData, fEventDoneData->fData, fEventDoneData->fDataSize );
1649       offset += fEventDoneData->fDataSize;
1650     }
1651     memcpy( reinterpret_cast<AliHLTUInt8_t*>(newEDD->fData)+offset, edd.fData, edd.fDataSize );
1652     if (fEventDoneData)
1653       delete [] reinterpret_cast<AliHLTUInt8_t*>( fEventDoneData );
1654     fEventDoneData = newEDD;
1655     fEventDoneDataSize = newSize;
1656   }
1657   else {
1658     memcpy( reinterpret_cast<AliHLTUInt8_t*>(fEventDoneData->fData)+fEventDoneData->fDataSize, edd.fData, edd.fDataSize );
1659     fEventDoneData->fDataSize += edd.fDataSize;
1660   }
1661   return iResult;
1662 }
1663
1664 int AliHLTComponent::ProcessEvent( const AliHLTComponentEventData& evtData,
1665                                    const AliHLTComponentBlockData* blocks, 
1666                                    AliHLTComponentTriggerData& trigData,
1667                                    AliHLTUInt8_t* outputPtr, 
1668                                    AliHLTUInt32_t& size,
1669                                    AliHLTUInt32_t& outputBlockCnt, 
1670                                    AliHLTComponentBlockData*& outputBlocks,
1671                                    AliHLTComponentEventDoneData*& edd )
1672 {
1673   // see header file for function documentation
1674   HLTLogKeyword(fChainId.c_str());
1675   ALIHLTCOMPONENT_BASE_STOPWATCH();
1676   int iResult=0;
1677   fCurrentEvent=evtData.fEventID;
1678   fCurrentEventData=evtData;
1679   fpInputBlocks=blocks;
1680   fCurrentInputBlock=-1;
1681   fSearchDataType=kAliHLTAnyDataType;
1682   fpOutputBuffer=outputPtr;
1683   fOutputBufferSize=size;
1684   fOutputBufferFilled=0;
1685   fOutputBlocks.clear();
1686   outputBlockCnt=0;
1687   outputBlocks=NULL;
1688
1689   AliHLTComponentBlockDataList forwardedBlocks;
1690
1691   // optional component statistics
1692   AliHLTComponentStatisticsList compStats;
1693   bool bAddComponentTableEntry=false;
1694   vector<AliHLTUInt32_t> parentComponentTables;
1695 #if defined(HLT_COMPONENT_STATISTICS)
1696   if ((fFlags&kDisableComponentStat)==0) {
1697     AliHLTComponentStatistics outputStat;
1698     memset(&outputStat, 0, sizeof(AliHLTComponentStatistics));
1699     outputStat.fStructSize=sizeof(AliHLTComponentStatistics);
1700     outputStat.fId=fChainIdCrc;
1701     if (fpBenchmark) {
1702       fpBenchmark->Stop();
1703       outputStat.fComponentCycleTime=(AliHLTUInt32_t)(fpBenchmark->RealTime()*ALIHLTCOMPONENT_STATTIME_SCALER);
1704       fpBenchmark->Reset();
1705       fpBenchmark->Start();
1706     }
1707     compStats.push_back(outputStat);
1708   }
1709 #endif // HLT_COMPONENT_STATISTICS
1710
1711   // data processing is skipped
1712   // -  if there are only steering events in the block list.
1713   //    For the sake of data source components data processing
1714   //    is not skipped if there is no block list at all or if it
1715   //    just contains the eventType block
1716   // - always skipped if the event is of type
1717   //   - gkAliEventTypeConfiguration
1718   //   - gkAliEventTypeReadPreprocessor
1719   const unsigned int skipModeDefault=0x1;
1720   const unsigned int skipModeForce=0x2;
1721   unsigned int bSkipDataProcessing=skipModeDefault;
1722
1723   // find special events
1724   if (fpInputBlocks && evtData.fBlockCnt>0) {
1725     // first look for all special events and execute in the appropriate
1726     // sequence afterwords
1727     int indexComConfEvent=-1;
1728     int indexUpdtDCSEvent=-1;
1729     int indexSOREvent=-1;
1730     int indexEOREvent=-1;
1731     int indexECSParamBlock=-1;
1732     for (unsigned int i=0; i<evtData.fBlockCnt && iResult>=0; i++) {
1733       if (fpInputBlocks[i].fDataType==kAliHLTDataTypeSOR) {
1734         indexSOREvent=i;
1735         // the AliHLTCalibrationProcessor relies on the SOR and EOR events
1736         bSkipDataProcessing&=~skipModeDefault;
1737       } else if (fpInputBlocks[i].fDataType==kAliHLTDataTypeRunType) {
1738         // run type string
1739         // handling is not clear yet
1740         if (fpInputBlocks[i].fPtr) {
1741           HLTDebug("got run type \"%s\"\n", fpInputBlocks[i].fPtr);
1742         }
1743       } else if (fpInputBlocks[i].fDataType==kAliHLTDataTypeEOR) {
1744         indexEOREvent=i;
1745         // the calibration processor relies on the SOR and EOR events
1746         bSkipDataProcessing&=~skipModeDefault;
1747       } else if (fpInputBlocks[i].fDataType==kAliHLTDataTypeDDL) {
1748         // DDL list
1749         // this event is most likely deprecated
1750       } else if (fpInputBlocks[i].fDataType==kAliHLTDataTypeComConf) {
1751         indexComConfEvent=i;
1752       } else if (fpInputBlocks[i].fDataType==kAliHLTDataTypeUpdtDCS) {
1753         indexUpdtDCSEvent=i;
1754       } else if (fpInputBlocks[i].fDataType==kAliHLTDataTypeEvent) {
1755         fEventType=fpInputBlocks[i].fSpecification;
1756
1757         // skip always in case of gkAliEventTypeConfiguration
1758         if (fpInputBlocks[i].fSpecification==gkAliEventTypeConfiguration) bSkipDataProcessing|=skipModeForce;
1759
1760         // skip always in case of gkAliEventTypeReadPreprocessor
1761         if (fpInputBlocks[i].fSpecification==gkAliEventTypeReadPreprocessor) bSkipDataProcessing|=skipModeForce;
1762
1763         // never skip if the event type block is the only block
1764         if (evtData.fBlockCnt==1) bSkipDataProcessing&=~skipModeDefault;
1765
1766       } else if (fpInputBlocks[i].fDataType==kAliHLTDataTypeComponentStatistics) {
1767         if (compStats.size()>0) {
1768           AliHLTUInt8_t* pData=reinterpret_cast<AliHLTUInt8_t*>(fpInputBlocks[i].fPtr);
1769           for (AliHLTUInt32_t offset=0;
1770                offset+sizeof(AliHLTComponentStatistics)<=fpInputBlocks[i].fSize;
1771                offset+=sizeof(AliHLTComponentStatistics)) {
1772             AliHLTComponentStatistics* pStat=reinterpret_cast<AliHLTComponentStatistics*>(pData+offset);
1773             if (pStat && compStats[0].fLevel<=pStat->fLevel) {
1774               compStats[0].fLevel=pStat->fLevel+1;
1775             }
1776             compStats.push_back(*pStat);
1777           }
1778         }
1779       } else if (fpInputBlocks[i].fDataType==kAliHLTDataTypeComponentTable) {
1780         forwardedBlocks.push_back(fpInputBlocks[i]);
1781         parentComponentTables.push_back(fpInputBlocks[i].fSpecification);
1782       } else if (fpInputBlocks[i].fDataType==kAliHLTDataTypeECSParam) {
1783         indexECSParamBlock=i;
1784       } else {
1785         // the processing function is called if there is at least one
1786         // non-steering data block. Steering blocks are not filtered out
1787         // for sake of performance 
1788         bSkipDataProcessing&=~skipModeDefault;
1789         if (compStats.size()>0) {
1790           compStats[0].fInputBlockCount++;
1791           compStats[0].fTotalInputSize+=fpInputBlocks[i].fSize;
1792         }
1793       }
1794     }
1795
1796     if (indexSOREvent>=0) {
1797       // start of run
1798       bAddComponentTableEntry=true;
1799       if (fpRunDesc==NULL) {
1800         fpRunDesc=new AliHLTRunDesc;
1801         if (fpRunDesc) *fpRunDesc=kAliHLTVoidRunDesc;
1802       }
1803       if (fpRunDesc) {
1804         AliHLTRunDesc rundesc;
1805         if ((iResult=CopyStruct(&rundesc, sizeof(AliHLTRunDesc), indexSOREvent, "AliHLTRunDesc", "SOR"))>0) {
1806           if (fpRunDesc->fRunNo==kAliHLTVoidRunNo) {
1807             *fpRunDesc=rundesc;
1808             HLTDebug("set run decriptor, run no %d", fpRunDesc->fRunNo);
1809             SetCDBRunNo(fpRunDesc->fRunNo);
1810           } else if (fpRunDesc->fRunNo!=rundesc.fRunNo) {
1811             HLTWarning("already set run properties run no %d, ignoring SOR with run no %d", fpRunDesc->fRunNo, rundesc.fRunNo);
1812           }
1813         }
1814       } else {
1815         iResult=-ENOMEM;
1816       }
1817
1818       if (indexECSParamBlock>=0) {
1819         if (fpInputBlocks[indexECSParamBlock].fSize>0) {
1820           const char* param=reinterpret_cast<const char*>(fpInputBlocks[indexECSParamBlock].fPtr);
1821           TString paramString;
1822           if (param[fpInputBlocks[indexECSParamBlock].fSize-1]!=0) {
1823             HLTWarning("ECS parameter string not terminated");
1824             paramString.Insert(0, param, fpInputBlocks[indexECSParamBlock].fSize);
1825             paramString+="";
1826           } else {
1827             paramString=param;
1828           }
1829           ScanECSParam(paramString.Data());
1830         } else {
1831           HLTWarning("empty ECS parameter received");
1832         }
1833       } else {
1834         // TODO: later on we might throw a warning here since the CTP trigger classes
1835         // should be mandatory
1836       }
1837     }
1838     if (indexEOREvent>=0) {
1839       fLastPushBackTime=0; // always send at EOR
1840       bAddComponentTableEntry=true;
1841       if (fpRunDesc!=NULL) {
1842         if (fpRunDesc) {
1843           AliHLTRunDesc rundesc;
1844           if ((iResult=CopyStruct(&rundesc, sizeof(AliHLTRunDesc), indexEOREvent, "AliHLTRunDesc", "SOR"))>0) {
1845             if (fpRunDesc->fRunNo!=rundesc.fRunNo) {
1846               HLTWarning("run no mismatch: SOR %d, EOR %d", fpRunDesc->fRunNo, rundesc.fRunNo);
1847             } else {
1848               HLTDebug("EOR run no %d", fpRunDesc->fRunNo);
1849             }
1850           }
1851           // we do not unload the fpRunDesc struct here in order to have the run information
1852           // available during the event processing
1853           // https://savannah.cern.ch/bugs/?39711
1854           // the info will be cleared in DeInit
1855         }
1856       } else {
1857         HLTWarning("did not receive SOR, ignoring EOR");
1858       }
1859     }
1860     if (indexComConfEvent>=0 || fEventType==gkAliEventTypeConfiguration) {
1861       TString cdbEntry;
1862       if (indexComConfEvent>=0 && fpInputBlocks[indexComConfEvent].fPtr!=NULL && fpInputBlocks[indexComConfEvent].fSize>0) {
1863         cdbEntry.Append(reinterpret_cast<const char*>(fpInputBlocks[indexComConfEvent].fPtr), fpInputBlocks[indexComConfEvent].fSize);
1864       }
1865       HLTDebug("received component configuration command: entry %s", cdbEntry.IsNull()?"none":cdbEntry.Data());
1866       int tmpResult=Reconfigure(cdbEntry[0]==0?NULL:cdbEntry.Data(), fChainId.c_str());
1867       if (tmpResult<0) {
1868         HLTWarning("reconfiguration of component %p (%s) failed with error code %d", this, GetComponentID(), tmpResult);
1869       }
1870     }
1871     if (indexUpdtDCSEvent>=0 || fEventType==gkAliEventTypeReadPreprocessor) {
1872       TString modules;
1873       if (fpInputBlocks[indexUpdtDCSEvent].fPtr!=NULL && fpInputBlocks[indexUpdtDCSEvent].fSize>0) {
1874         modules.Append(reinterpret_cast<const char*>(fpInputBlocks[indexUpdtDCSEvent].fPtr), fpInputBlocks[indexUpdtDCSEvent].fSize);
1875       }
1876       HLTDebug("received preprocessor update command: detectors %s", modules.IsNull()?"ALL":modules.Data());
1877       int tmpResult=ReadPreprocessorValues(modules[0]==0?"ALL":modules.Data());
1878       if (tmpResult<0) {
1879         HLTWarning("preprocessor update of component %p (%s) failed with error code %d", this, GetComponentID(), tmpResult);
1880       }
1881     }
1882   } else {
1883     // processing function needs to be called if there are no input data
1884     // blocks in order to make data source components working.
1885     bSkipDataProcessing&=~skipModeDefault;
1886   }
1887
1888   // data processing is not skipped if the component explicitly asks
1889   // for the private blocks
1890   if ((fFlags&kRequireSteeringBlocks)!=0) bSkipDataProcessing=0;
1891
1892   if (fpCTPData) {
1893     // set the active triggers for this event
1894     fpCTPData->SetTriggers(trigData);
1895     // increment CTP trigger counters if available
1896     if (IsDataEvent()) fpCTPData->Increment(trigData);
1897   }
1898
1899   AliHLTComponentBlockDataList blockData;
1900   if (iResult>=0 && !bSkipDataProcessing)
1901   { // dont delete, sets the scope for the stopwatch guard
1902     // do not use ALIHLTCOMPONENT_DA_STOPWATCH(); macro
1903     // in order to avoid 'shadowed variable' warning
1904     AliHLTStopwatchGuard swguard2(fpStopwatches!=NULL?reinterpret_cast<TStopwatch*>(fpStopwatches->At((int)kSWDA)):NULL);
1905     iResult=DoProcessing(evtData, blocks, trigData, outputPtr, size, blockData, edd);
1906   } // end of the scope of the stopwatch guard
1907   if (iResult>=0 && !bSkipDataProcessing) {
1908     if (fOutputBlocks.size()>0) {
1909       // High Level interface
1910
1911       //HLTDebug("got %d block(s) via high level interface", fOutputBlocks.size());      
1912       // sync memory files and descriptors
1913       AliHLTMemoryFilePList::iterator element=fMemFiles.begin();
1914       int i=0;
1915       while (element!=fMemFiles.end() && iResult>=0) {
1916         if (*element) {
1917           if ((*element)->IsClosed()==0) {
1918             HLTWarning("memory file has not been closed, force flush");
1919             iResult=CloseMemoryFile(*element);
1920           }
1921         }
1922         element++; i++;
1923       }
1924
1925       if (iResult>=0) {
1926         // create the descriptor list
1927         if (blockData.size()>0) {
1928           HLTError("low level and high interface must not be mixed; use PushBack methods to insert data blocks");
1929           iResult=-EFAULT;
1930         } else {
1931           if (compStats.size()>0 && IsDataEvent()) {
1932             int offset=AddComponentStatistics(fOutputBlocks, fpOutputBuffer, fOutputBufferSize, fOutputBufferFilled, compStats);
1933             if (offset>0) fOutputBufferFilled+=offset;
1934           }
1935           if (bAddComponentTableEntry) {
1936             int offset=AddComponentTableEntry(fOutputBlocks, fpOutputBuffer, fOutputBufferSize, fOutputBufferFilled, parentComponentTables);
1937             if (offset>0) size+=offset;
1938           }
1939           if (forwardedBlocks.size()>0) {
1940             fOutputBlocks.insert(fOutputBlocks.end(), forwardedBlocks.begin(), forwardedBlocks.end());
1941           }
1942           iResult=MakeOutputDataBlockList(fOutputBlocks, &outputBlockCnt, &outputBlocks);
1943           size=fOutputBufferFilled;
1944         }
1945       }
1946     } else {
1947       // Low Level interface
1948       if (compStats.size()>0) {
1949         int offset=AddComponentStatistics(blockData, fpOutputBuffer, fOutputBufferSize, size, compStats);
1950         if (offset>0) size+=offset;
1951       }
1952       if (bAddComponentTableEntry) {
1953         int offset=AddComponentTableEntry(blockData, fpOutputBuffer, fOutputBufferSize, size, parentComponentTables);
1954         if (offset>0) size+=offset;
1955       }
1956       if (forwardedBlocks.size()>0) {
1957         blockData.insert(blockData.end(), forwardedBlocks.begin(), forwardedBlocks.end());
1958       }
1959       iResult=MakeOutputDataBlockList(blockData, &outputBlockCnt, &outputBlocks);
1960     }
1961     if (iResult<0) {
1962       HLTFatal("component %s (%p): can not convert output block descriptor list", GetComponentID(), this);
1963     }
1964   }
1965   if (iResult<0 || bSkipDataProcessing) {
1966     outputBlockCnt=0;
1967     outputBlocks=NULL;
1968   }
1969   CleanupInputObjects();
1970   if (iResult>=0 && IsDataEvent()) {
1971     IncrementEventCounter();
1972   }
1973   if (outputBlockCnt==0) {
1974     // no output blocks, set size to 0
1975     size=0;
1976   }
1977
1978   // reset the internal EventData struct
1979   FillEventData(fCurrentEventData);
1980
1981   // reset the active triggers
1982   if (fpCTPData) fpCTPData->SetTriggers(0);
1983
1984   // set the time for the pushback period
1985   if (fPushbackPeriod>0) {
1986     // suppress the output
1987     TDatime time;
1988     if (fLastPushBackTime<0 || (int)time.Get()-fLastPushBackTime>=fPushbackPeriod) {
1989       fLastPushBackTime=time.Get();
1990     }
1991   }
1992
1993   return iResult;
1994 }
1995
1996 int  AliHLTComponent::AddComponentStatistics(AliHLTComponentBlockDataList& blocks, 
1997                                              AliHLTUInt8_t* buffer,
1998                                              AliHLTUInt32_t bufferSize,
1999                                              AliHLTUInt32_t offset,
2000                                              AliHLTComponentStatisticsList& stats) const
2001 {
2002   // see header file for function documentation
2003   int iResult=0;
2004   if ((fFlags&kDisableComponentStat)!=0) return 0;
2005 #if defined(HLT_COMPONENT_STATISTICS)
2006   if (stats.size()==0) return -ENOENT;
2007   // check if there is space for at least one entry
2008   if (offset+sizeof(AliHLTComponentStatistics)>bufferSize) return 0;
2009   stats[0].fTotalOutputSize=offset;
2010   stats[0].fOutputBlockCount=blocks.size();
2011   if (fpBenchmark) {
2012     fpBenchmark->Stop();
2013     stats[0].fTime=(AliHLTUInt32_t)(fpBenchmark->RealTime()*ALIHLTCOMPONENT_STATTIME_SCALER);
2014     stats[0].fCTime=(AliHLTUInt32_t)(fpBenchmark->CpuTime()*ALIHLTCOMPONENT_STATTIME_SCALER);
2015     fpBenchmark->Continue();
2016   }
2017   if (offset+stats.size()*sizeof(AliHLTComponentStatistics)>bufferSize) {
2018     AliHLTUInt32_t removedLevel=0;
2019     do {
2020       // remove all entries of the level of the last entry
2021       removedLevel=stats.back().fLevel;
2022       AliHLTComponentStatisticsList::iterator element=stats.begin();
2023       element++;
2024       while (element!=stats.end()) {
2025         if (element->fLevel<=removedLevel) {
2026           element=stats.erase(element);
2027         } else {
2028           element++;
2029         }
2030       }
2031     } while (stats.size()>1 && 
2032              (offset+stats.size()*sizeof(AliHLTComponentStatistics)>bufferSize));
2033   }
2034   assert(stats.size()>0);
2035   if (stats.size()==0) return 0;
2036
2037   if (offset+stats.size()*sizeof(AliHLTComponentStatistics)<=bufferSize) {
2038     AliHLTComponentBlockData bd;
2039     FillBlockData( bd );
2040     bd.fOffset        = offset;
2041     bd.fSize          = stats.size()*sizeof(AliHLTComponentStatistics);
2042     bd.fDataType      = kAliHLTDataTypeComponentStatistics;
2043     bd.fSpecification = kAliHLTVoidDataSpec;
2044     unsigned int master=0;
2045     for (unsigned int i=1; i<blocks.size(); i++) {
2046       if (blocks[i].fSize>blocks[master].fSize && 
2047           !MatchExactly(blocks[i].fDataType, kAliHLTVoidDataType|kAliHLTDataOriginPrivate))
2048         master=i;
2049     }
2050     if (blocks.size()>0 && !MatchExactly(blocks[master].fDataType, kAliHLTVoidDataType|kAliHLTDataOriginPrivate)) {
2051       // take the data origin of the biggest block as specification
2052       // this is similar to the treatment in the HOMER interface. For traditional
2053       // reasons, the bytes are swapped there on a little endian architecture, so
2054       // we do it as well.
2055       memcpy(&bd.fSpecification, &blocks[master].fDataType.fOrigin, sizeof(bd.fSpecification));
2056 #ifdef R__BYTESWAP // set on little endian architectures
2057       bd.fSpecification=((bd.fSpecification & 0xFFULL) << 24) | 
2058         ((bd.fSpecification & 0xFF00ULL) << 8) | 
2059         ((bd.fSpecification & 0xFF0000ULL) >> 8) | 
2060         ((bd.fSpecification & 0xFF000000ULL) >> 24);
2061 #endif
2062     }
2063     memcpy(buffer+offset, &(stats[0]), bd.fSize);
2064     blocks.push_back(bd);
2065     iResult=bd.fSize;
2066   }
2067 #else
2068   if (blocks.size() && buffer && bufferSize && offset && stats.size()) {
2069     // get rid of warning
2070   }
2071 #endif // HLT_COMPONENT_STATISTICS
2072   return iResult;
2073 }
2074
2075 int  AliHLTComponent::AddComponentTableEntry(AliHLTComponentBlockDataList& blocks, 
2076                                              AliHLTUInt8_t* buffer,
2077                                              AliHLTUInt32_t bufferSize,
2078                                              AliHLTUInt32_t offset,
2079                                              const vector<AliHLTUInt32_t>& parents) const
2080 {
2081   // see header file for function documentation
2082   int iResult=0;
2083   if ((fFlags&kDisableComponentStat)!=0) return 0;
2084 #if defined(HLT_COMPONENT_STATISTICS)
2085   // the payload consists of the AliHLTComponentTableEntry struct,
2086   // followed by a an array of 32bit crc chain ids and the component
2087   // description string
2088   unsigned int payloadSize=sizeof(AliHLTComponentTableEntry);
2089   payloadSize+=parents.size()*sizeof(AliHLTUInt32_t);
2090
2091   // the component description has the following format:
2092   // chain-id{component-id:arguments}
2093   const char* componentId=const_cast<AliHLTComponent*>(this)->GetComponentID();
2094   unsigned int descriptionSize=fChainId.size()+1;
2095   descriptionSize+=2; // the '{}' around the component id
2096   descriptionSize+=strlen(componentId);
2097   descriptionSize+=1; // the ':' between component id and arguments
2098   descriptionSize+=fComponentArgs.size();
2099
2100   payloadSize+=descriptionSize;
2101   if (buffer && (offset+payloadSize<=bufferSize)) {
2102     AliHLTUInt8_t* pTgt=buffer+offset;
2103     memset(pTgt, 0, payloadSize);
2104
2105     // write entry
2106     AliHLTComponentTableEntry* pEntry=reinterpret_cast<AliHLTComponentTableEntry*>(pTgt);
2107     pEntry->fStructSize=sizeof(AliHLTComponentTableEntry);
2108     pEntry->fNofParents=parents.size();
2109     pEntry->fSizeDescription=descriptionSize;
2110     pTgt=pEntry->fBuffer;
2111
2112     // write array of parents
2113     if (parents.size()>0) {
2114       unsigned int copy=parents.size()*sizeof(vector<AliHLTUInt32_t>::value_type);
2115       memcpy(pTgt, &parents[0], parents.size()*sizeof(vector<AliHLTUInt32_t>::value_type));
2116       pTgt+=copy;
2117     }
2118
2119     // write component description
2120     memcpy(pTgt, fChainId.c_str(), fChainId.size());
2121     pTgt+=fChainId.size();
2122     *pTgt++='{';
2123     memcpy(pTgt, componentId, strlen(componentId));
2124     pTgt+=strlen(componentId);
2125     *pTgt++=':';
2126     memcpy(pTgt, fComponentArgs.c_str(), fComponentArgs.size());
2127     pTgt+=fComponentArgs.size();
2128     *pTgt++='}';
2129     *pTgt++=0;
2130
2131     AliHLTComponentBlockData bd;
2132     FillBlockData( bd );
2133     bd.fOffset        = offset;
2134     bd.fSize          = payloadSize;
2135     bd.fDataType      = kAliHLTDataTypeComponentTable;
2136     bd.fSpecification = fChainIdCrc;
2137     blocks.push_back(bd);
2138     iResult=bd.fSize;
2139   }
2140 #else
2141   if (blocks.size() && buffer && bufferSize && offset && parents.size()) {
2142     // get rid of warning
2143   }
2144 #endif // HLT_COMPONENT_STATISTICS
2145   return iResult;
2146 }
2147
2148 AliHLTComponent::AliHLTStopwatchGuard::AliHLTStopwatchGuard()
2149   :
2150   fpStopwatch(NULL),
2151   fpPrec(NULL)
2152 {
2153   // standard constructor (not for use)
2154 }
2155
2156 AliHLTComponent::AliHLTStopwatchGuard::AliHLTStopwatchGuard(TStopwatch* pStopwatch)
2157   :
2158   fpStopwatch(pStopwatch),
2159   fpPrec(NULL)
2160 {
2161   // constructor
2162
2163   // check for already existing guard
2164   if (fgpCurrent) fpPrec=fgpCurrent;
2165   fgpCurrent=this;
2166
2167   // stop the preceeding guard if it controls a different stopwatch
2168   int bStart=1;
2169   if (fpPrec && fpPrec!=this) bStart=fpPrec->Hold(fpStopwatch);
2170
2171   // start the stopwatch if the current guard controls a different one
2172   if (fpStopwatch && bStart==1) fpStopwatch->Start(kFALSE);
2173 }
2174
2175 AliHLTComponent::AliHLTStopwatchGuard::AliHLTStopwatchGuard(const AliHLTStopwatchGuard&)
2176   :
2177   fpStopwatch(NULL),
2178   fpPrec(NULL)
2179 {
2180   //
2181   // copy constructor not for use
2182   //
2183 }
2184
2185 AliHLTComponent::AliHLTStopwatchGuard& AliHLTComponent::AliHLTStopwatchGuard::operator=(const AliHLTStopwatchGuard&)
2186 {
2187   //
2188   // assignment operator not for use
2189   //
2190   fpStopwatch=NULL;
2191   fpPrec=NULL;
2192   return *this;
2193 }
2194
2195 AliHLTComponent::AliHLTStopwatchGuard* AliHLTComponent::AliHLTStopwatchGuard::fgpCurrent=NULL;
2196
2197 AliHLTComponent::AliHLTStopwatchGuard::~AliHLTStopwatchGuard()
2198 {
2199   // destructor
2200
2201   // resume the preceeding guard if it controls a different stopwatch
2202   int bStop=1;
2203   if (fpPrec && fpPrec!=this) bStop=fpPrec->Resume(fpStopwatch);
2204
2205   // stop the stopwatch if the current guard controls a different one
2206   if (fpStopwatch && bStop==1) fpStopwatch->Stop();
2207
2208   // resume to the preceeding guard
2209   fgpCurrent=fpPrec;
2210 }
2211
2212 int AliHLTComponent::AliHLTStopwatchGuard::Hold(const TStopwatch* pSucc)
2213 {
2214   // see header file for function documentation
2215   if (fpStopwatch!=NULL && fpStopwatch!=pSucc) fpStopwatch->Stop();
2216   return fpStopwatch!=pSucc?1:0;
2217 }
2218
2219 int AliHLTComponent::AliHLTStopwatchGuard::Resume(const TStopwatch* pSucc)
2220 {
2221   // see header file for function documentation
2222   if (fpStopwatch!=NULL && fpStopwatch!=pSucc) fpStopwatch->Start(kFALSE);
2223   return fpStopwatch!=pSucc?1:0;
2224 }
2225
2226 int AliHLTComponent::SetStopwatch(TObject* pSW, AliHLTStopwatchType type) 
2227 {
2228   // see header file for function documentation
2229   int iResult=0;
2230   if (pSW!=NULL && type<kSWTypeCount) {
2231     if (fpStopwatches) {
2232       TObject* pObj=fpStopwatches->At((int)type);
2233       if (pSW==NULL        // explicit reset
2234           || pObj==NULL) { // explicit set
2235         fpStopwatches->AddAt(pSW, (int)type);
2236       } else if (pObj!=pSW) {
2237         HLTWarning("stopwatch %d already set, reset first", (int)type);
2238         iResult=-EBUSY;
2239       }
2240     }
2241   } else {
2242     iResult=-EINVAL;
2243   }
2244   return iResult;
2245 }
2246
2247 int AliHLTComponent::SetStopwatches(TObjArray* pStopwatches)
2248 {
2249   // see header file for function documentation
2250   if (pStopwatches==NULL) return -EINVAL;
2251
2252   int iResult=0;
2253   for (int i=0 ; i<(int)kSWTypeCount && pStopwatches->GetEntriesFast(); i++)
2254     SetStopwatch(pStopwatches->At(i), (AliHLTStopwatchType)i);
2255   return iResult;
2256 }
2257
2258 AliHLTUInt32_t AliHLTComponent::GetRunNo() const
2259 {
2260   // see header file for function documentation
2261
2262   // 2010-02-11 OCDB is now the reliable source for the run number, it is
2263   // initialized either by aliroot or the external interface depending
2264   // on the environment. It turned out that the rundescriptor is not set
2265   // in the aliroot mode, resulting in an invalid run number. However this
2266   // did not cause problems until now. OCDB initialization has been revised
2267   // already in 10/2009. This was a remnant.
2268   // Have to check whether we get rid of the rundescriptor at some point.
2269   if (fpRunDesc==NULL) return AliHLTMisc::Instance().GetCDBRunNo();
2270   if (fpRunDesc->fRunNo!=(unsigned)AliHLTMisc::Instance().GetCDBRunNo()) {
2271     HLTWarning("run number mismatch: ocdb %d   run descriptor %d", AliHLTMisc::Instance().GetCDBRunNo(), fpRunDesc->fRunNo);
2272   }
2273   return fpRunDesc->fRunNo;
2274 }
2275
2276 AliHLTUInt32_t AliHLTComponent::GetRunType() const
2277 {
2278   // see header file for function documentation
2279   if (fpRunDesc==NULL) return kAliHLTVoidRunType;
2280   return fpRunDesc->fRunType;
2281 }
2282
2283 AliHLTUInt32_t    AliHLTComponent::GetTimeStamp() const
2284 {
2285   // see header file for function documentation
2286   if (fCurrentEventData.fEventCreation_s) {
2287     return  fCurrentEventData.fEventCreation_s;
2288   }
2289   // using the actual UTC if the time stamp was not set by the framework
2290   return static_cast<AliHLTUInt32_t>(time(NULL));
2291 }
2292
2293 AliHLTUInt32_t    AliHLTComponent::GetPeriodNumber() const
2294 {
2295   // see header file for function documentation
2296   return (GetEventId()>>36)&0xfffffff;
2297 }
2298
2299 AliHLTUInt32_t    AliHLTComponent::GetOrbitNumber() const
2300 {
2301   // see header file for function documentation
2302   return (GetEventId()>>12)&0xffffff;
2303 }
2304
2305 AliHLTUInt16_t    AliHLTComponent::GetBunchCrossNumber() const
2306 {
2307   // see header file for function documentation
2308   return GetEventId()&0xfff;
2309 }
2310
2311 bool AliHLTComponent::IsDataEvent(AliHLTUInt32_t* pTgt) const
2312 {
2313   // see header file for function documentation
2314   if (pTgt) *pTgt=fEventType;
2315   return (fEventType==gkAliEventTypeData ||
2316           fEventType==gkAliEventTypeDataReplay ||
2317           fEventType==gkAliEventTypeCalibration);
2318 }
2319
2320 int AliHLTComponent::CopyStruct(void* pStruct, unsigned int iStructSize, unsigned int iBlockNo,
2321                                 const char* structname, const char* eventname)
2322 {
2323   // see header file for function documentation
2324   int iResult=0;
2325   if (pStruct!=NULL && iStructSize>sizeof(AliHLTUInt32_t)) {
2326     if (fpInputBlocks!=NULL && iBlockNo<fCurrentEventData.fBlockCnt) {
2327       AliHLTUInt32_t* pTgt=(AliHLTUInt32_t*)pStruct;
2328       if (fpInputBlocks[iBlockNo].fPtr && fpInputBlocks[iBlockNo].fSize) {
2329         AliHLTUInt32_t copy=*((AliHLTUInt32_t*)fpInputBlocks[iBlockNo].fPtr);
2330         if (fpInputBlocks[iBlockNo].fSize!=copy) {
2331           HLTWarning("%s event: mismatch of block size (%d) and structure size (%d)", eventname, fpInputBlocks[iBlockNo].fSize, copy);
2332           if (copy>fpInputBlocks[iBlockNo].fSize) copy=fpInputBlocks[iBlockNo].fSize;
2333         }
2334         if (copy!=iStructSize) {
2335           HLTWarning("%s event: mismatch in %s version (data type version %d)", eventname, structname, ALIHLT_DATA_TYPES_VERSION);
2336           if (copy>iStructSize) {
2337             copy=iStructSize;
2338           } else {
2339             memset(pTgt, 0, iStructSize);
2340           }
2341         }
2342         memcpy(pTgt, fpInputBlocks[iBlockNo].fPtr, copy);
2343         *pTgt=iStructSize;
2344         iResult=copy;
2345       } else {
2346         HLTWarning("%s event: missing data block", eventname);
2347       }
2348     } else {
2349       iResult=-ENODATA;
2350     }
2351   } else {
2352     HLTError("invalid struct");
2353     iResult=-EINVAL;
2354   }
2355   return iResult;
2356 }
2357
2358 AliHLTUInt32_t AliHLTComponent::CalculateChecksum(const AliHLTUInt8_t* buffer, int size)
2359 {
2360   // see header file for function documentation
2361   AliHLTUInt32_t  remainder = 0; 
2362   const AliHLTUInt8_t crcwidth=(8*sizeof(AliHLTUInt32_t));
2363   const AliHLTUInt32_t topbit=1 << (crcwidth-1);
2364   const AliHLTUInt32_t polynomial=0xD8;  /* 11011 followed by 0's */
2365
2366   // code from
2367   // http://www.netrino.com/Embedded-Systems/How-To/CRC-Calculation-C-Code
2368
2369   /*
2370    * Perform modulo-2 division, a byte at a time.
2371    */
2372   for (int byte = 0; byte < size; ++byte)
2373     {
2374       /*
2375        * Bring the next byte into the remainder.
2376        */
2377       remainder ^= (buffer[byte] << (crcwidth - 8));
2378
2379       /*
2380        * Perform modulo-2 division, a bit at a time.
2381        */
2382       for (uint8_t bit = 8; bit > 0; --bit)
2383         {
2384           /*
2385            * Try to divide the current data bit.
2386            */
2387           if (remainder & topbit)
2388             {
2389               remainder = (remainder << 1) ^ polynomial;
2390             }
2391           else
2392             {
2393               remainder = (remainder << 1);
2394             }
2395         }
2396     }
2397
2398   /*
2399    * The final remainder is the CRC result.
2400    */
2401   return (remainder);
2402 }
2403
2404 int AliHLTComponent::ExtractComponentTableEntry(const AliHLTUInt8_t* pBuffer, AliHLTUInt32_t size,
2405                                                 string& retChainId, string& retCompId, string& retCompArgs,
2406                                                 vector<AliHLTUInt32_t>& parents)
2407 {
2408   // see header file for function documentation
2409   retChainId.clear();
2410   retCompId.clear();
2411   retCompArgs.clear();
2412   parents.clear();
2413   if (!pBuffer || size==0) return 0;
2414
2415   const AliHLTComponentTableEntry* pEntry=reinterpret_cast<const AliHLTComponentTableEntry*>(pBuffer);
2416   if (size<8/* the initial size of the structure*/ ||
2417       pEntry==NULL || pEntry->fStructSize<8) return -ENOMSG;
2418   const AliHLTUInt32_t* pParents=reinterpret_cast<const AliHLTUInt32_t*>(pEntry->fBuffer);
2419   const AliHLTUInt8_t* pEnd=pBuffer+size;
2420
2421   if (pParents+pEntry->fNofParents>=reinterpret_cast<const AliHLTUInt32_t*>(pEnd)) return -ENODEV;
2422   for (unsigned int i=0; i<pEntry->fNofParents; i++, pParents++) {
2423     parents.push_back(*pParents);
2424   }
2425
2426   const char* pDescription=reinterpret_cast<const char*>(pParents);
2427   if (pDescription+pEntry->fSizeDescription>=reinterpret_cast<const char*>(pEnd) ||
2428       *(pDescription+pEntry->fSizeDescription)!=0) {
2429     return -EBADF;
2430   }
2431
2432   TString descriptor=reinterpret_cast<const char*>(pDescription);
2433   TString chainId;
2434   TString compId;
2435   TString compArgs;
2436   TObjArray* pTokens=descriptor.Tokenize("{");
2437   if (pTokens) {
2438     int n=0;
2439     if (pTokens->GetEntriesFast()>n) {
2440       retChainId=((TObjString*)pTokens->At(n++))->GetString();
2441     }
2442     if (pTokens->GetEntriesFast()>n) {
2443       compId=((TObjString*)pTokens->At(n++))->GetString();
2444     }
2445     delete pTokens;
2446   }
2447   if (!compId.IsNull() && (pTokens=compId.Tokenize(":"))!=NULL) {
2448     int n=0;
2449     if (pTokens->GetEntriesFast()>n) {
2450       compId=((TObjString*)pTokens->At(n++))->GetString();
2451     }
2452     if (pTokens->GetEntriesFast()>n) {
2453       compArgs=((TObjString*)pTokens->At(n++))->GetString();
2454     }
2455     delete pTokens;
2456   }
2457   compId.ReplaceAll("}", "");
2458   compArgs.ReplaceAll("}", "");
2459
2460   retCompId=compId;
2461   retCompArgs=compArgs;
2462
2463   if (retChainId.size()==0) return -ENODATA;
2464
2465   return 1;
2466 }
2467
2468 int AliHLTComponent::ExtractTriggerData(
2469     const AliHLTComponentTriggerData& trigData,
2470     const AliHLTUInt8_t (**attributes)[gkAliHLTBlockDAttributeCount],
2471     AliHLTUInt64_t* status,
2472     const AliRawDataHeader** cdh,
2473     AliHLTReadoutList* readoutlist,
2474     bool printErrors
2475   )
2476 {
2477   // see header file for function documentation
2478   
2479   // Check that the trigger data structure is the correct size.
2480   if (trigData.fStructSize != sizeof(AliHLTComponentTriggerData))
2481   {
2482     if (printErrors)
2483     {
2484       AliHLTLogging log;
2485       log.LoggingVarargs(kHLTLogError, Class_Name(), FUNCTIONNAME(), __FILE__, __LINE__,
2486           "Invalid trigger structure size: %d but expected %d.", trigData.fStructSize, sizeof(AliHLTComponentTriggerData)
2487         );
2488     }
2489     return -ENOENT;
2490   }
2491   
2492   // Get the size of the AliHLTEventTriggerData structure without the readout list part.
2493   // The way we do this here should also handle memory alignment correctly.
2494   AliHLTEventTriggerData* dummy = NULL;
2495   size_t sizeWithoutReadout = (char*)(&dummy->fReadoutList) - (char*)(dummy);
2496   
2497   // Check that the trigger data pointer points to data of a size we can handle.
2498   // Either it is the size of AliHLTEventTriggerData or the size of the old
2499   // version of AliHLTEventTriggerData using AliHLTEventDDLV0.
2500   if (trigData.fDataSize != sizeof(AliHLTEventTriggerData) and
2501       trigData.fDataSize != sizeWithoutReadout + sizeof(AliHLTEventDDLV0)
2502      )
2503   {
2504     if (printErrors)
2505     {
2506       AliHLTLogging log;
2507       log.LoggingVarargs(kHLTLogError, Class_Name(), FUNCTIONNAME(), __FILE__, __LINE__,
2508           "Invalid trigger data size: %d but expected %d.", trigData.fDataSize, sizeof(AliHLTEventTriggerData)
2509         );
2510     }
2511     return -EBADF;
2512   }
2513   
2514   AliHLTEventTriggerData* evtData = reinterpret_cast<AliHLTEventTriggerData*>(trigData.fData);
2515   assert(evtData != NULL);
2516   
2517   // Check that the CDH has 8 words.
2518   if (cdh != NULL and evtData->fCommonHeaderWordCnt != 8)
2519   {
2520     if (printErrors)
2521     {
2522       AliHLTLogging log;
2523       log.LoggingVarargs(kHLTLogError, Class_Name(), FUNCTIONNAME(), __FILE__, __LINE__,
2524           "Common Data Header (CDH) has wrong number of data words: %d but expected %d",
2525           evtData->fCommonHeaderWordCnt, sizeof(AliRawDataHeader)/sizeof(AliHLTUInt32_t)
2526         );
2527     }
2528     return -EBADMSG;
2529   }
2530   
2531   // Check that the readout list has the correct count of words. i.e. something we can handle,
2532   if (readoutlist != NULL and
2533       evtData->fReadoutList.fCount != (unsigned)gkAliHLTDDLListSizeV0 and
2534       evtData->fReadoutList.fCount != (unsigned)gkAliHLTDDLListSizeV1
2535      )
2536   {
2537     if (printErrors)
2538     {
2539       AliHLTLogging log;
2540       log.LoggingVarargs(kHLTLogError, Class_Name(), FUNCTIONNAME(), __FILE__, __LINE__,
2541           "Readout list structure has wrong number of data words: %d but expected %d",
2542           evtData->fReadoutList.fCount, gkAliHLTDDLListSize
2543         );
2544     }
2545     return -EPROTO;
2546   }
2547   
2548   if (attributes != NULL)
2549   {
2550     *attributes = &evtData->fAttributes;
2551   }
2552   if (status != NULL)
2553   {
2554     *status = evtData->fHLTStatus;
2555   }
2556   if (cdh != NULL)
2557   {
2558     const AliRawDataHeader* cdhptr = reinterpret_cast<const AliRawDataHeader*>(&evtData->fCommonHeader);
2559     *cdh = cdhptr;
2560   }
2561   if (readoutlist != NULL)
2562   {
2563     *readoutlist = AliHLTReadoutList(evtData->fReadoutList);
2564   }
2565   return 0;
2566 }
2567
2568 int AliHLTComponent::LoggingVarargs(AliHLTComponentLogSeverity severity, 
2569                                     const char* originClass, const char* originFunc,
2570                                     const char* file, int line, ... ) const
2571 {
2572   // see header file for function documentation
2573   int iResult=0;
2574
2575   va_list args;
2576   va_start(args, line);
2577
2578   // logging function needs to be const in order to be called from const member functions
2579   // without problems. But at this point we face the problem with virtual members which
2580   // are not necessarily const.
2581   AliHLTComponent* nonconst=const_cast<AliHLTComponent*>(this);
2582   AliHLTLogging::SetLogString(this, ", %p", "%s (%s_pfmt_): ", 
2583                               fChainId[0]!=0?fChainId.c_str():nonconst->GetComponentID(),
2584                               nonconst->GetComponentID());
2585   iResult=SendMessage(severity, originClass, originFunc, file, line, AliHLTLogging::BuildLogString(NULL, args, true /*append*/));
2586   va_end(args);
2587
2588   return iResult;
2589 }
2590
2591 TUUID AliHLTComponent::GenerateGUID()
2592 {
2593   // Generates a globally unique identifier.
2594
2595   // Start by creating a new UUID. We cannot use the one automatically generated
2596   // by ROOT because the algorithm used will not guarantee unique IDs when generating
2597   // these UUIDs at a high rate in parallel.
2598   TUUID uuid;
2599   // We then use the generated UUID to form part of the random number seeds which
2600   // will be used to generate a proper random UUID. For good measure we use a MD5
2601   // hash also. Note that we want to use the TUUID class because it will combine the
2602   // host address information into the UUID. Using gSystem->GetHostByName() apparently
2603   // can cause problems on Windows machines with a firewall, because it always tries
2604   // to contact a DNS. The TUUID class handles this case appropriately.
2605   union
2606   {
2607     UChar_t buf[16];
2608     UShort_t word[8];
2609     UInt_t dword[4];
2610   };
2611   uuid.GetUUID(buf);
2612   TMD5 md5;
2613   md5.Update(buf, sizeof(buf));
2614   TMD5 md52 = md5;
2615   md5.Final(buf);
2616   dword[0] += gSystem->GetUid();
2617   dword[1] += gSystem->GetGid();
2618   dword[2] += gSystem->GetPid();
2619   for (int i = 0; i < 4; ++i)
2620   {
2621     gRandom->SetSeed(dword[i]);
2622     dword[i] = gRandom->Integer(0xFFFFFFFF);
2623   }
2624   md52.Update(buf, sizeof(buf));
2625   md52.Final(buf);
2626   // To keep to the standard we need to set the version and reserved bits.
2627   word[3] = (word[3] & 0x0FFF) | 0x4000;
2628   buf[8] = (buf[8] & 0x3F) | 0x80;
2629
2630   // Create the name of the new class and file.
2631   char uuidstr[64];
2632   sprintf(uuidstr, "%08x-%04x-%04x-%02x%02x-%02x%02x%02x%02x%02x%02x",
2633     dword[0], word[2], word[3], buf[8], buf[9], buf[10], buf[11], buf[12], buf[13], buf[14], buf[15]
2634   );
2635
2636   uuid.SetUUID(uuidstr);
2637   return uuid;
2638 }
2639
2640 int AliHLTComponent::ScanECSParam(const char* ecsParam)
2641 {
2642   // see header file for function documentation
2643
2644   // format of the parameter string from ECS
2645   // <command>;<parameterkey>=<parametervalue>;<parameterkey>=<parametervalue>;...
2646   // search for a subset of the parameterkeys
2647   //   RUN_TYPE=
2648   //   RUN_NUMBER=
2649   //   HLT_IN_DDL_LIST=
2650   //   CTP_TRIGGER_CLASS=
2651   //   DATA_FORMAT_VERSION=
2652   //   BEAM_TYPE=
2653   //   HLT_OUT_DDL_LIST=
2654   //   HLT_TRIGGER_CODE=
2655   //   DETECTOR_LIST=
2656   //   HLT_MODE=
2657   // The command apears not to be sent by the online framework
2658   int iResult=0;
2659   TString string=ecsParam;
2660   TObjArray* parameter=string.Tokenize(";");
2661   if (parameter) {
2662     for (int i=0; i<parameter->GetEntriesFast(); i++) {
2663       TString entry=((TObjString*)parameter->At(i))->GetString();
2664       HLTDebug("scanning ECS entry: %s", entry.Data());
2665       TObjArray* entryParams=entry.Tokenize("=");
2666       if (entryParams) {
2667         if (entryParams->GetEntriesFast()>1) {
2668           if ((((TObjString*)entryParams->At(0))->GetString()).CompareTo("CTP_TRIGGER_CLASS")==0) {
2669             int result=InitCTPTriggerClasses((((TObjString*)entryParams->At(1))->GetString()).Data());
2670             if (iResult>=0 && result<0) iResult=result;
2671           } else {
2672             // TODO: scan the other parameters
2673             // e.g. consistency check of run number
2674           }
2675         }
2676         delete entryParams;
2677       }
2678     }
2679     delete parameter;
2680   }
2681
2682   return iResult;
2683 }
2684
2685 int AliHLTComponent::SetupCTPData()
2686 {
2687   // see header file for function documentation
2688   if (fpCTPData) delete fpCTPData;
2689   fpCTPData=new AliHLTCTPData;
2690   if (!fpCTPData) return -ENOMEM;
2691   return 0;
2692 }
2693
2694 int AliHLTComponent::InitCTPTriggerClasses(const char* ctpString)
2695 {
2696   // see header file for function documentation
2697   if (!fpCTPData) return 0; // silently accept as the component has to announce that it want's the CTP info
2698   return fpCTPData->InitCTPTriggerClasses(ctpString);
2699 }
2700
2701 bool AliHLTComponent::EvaluateCTPTriggerClass(const char* expression, AliHLTComponentTriggerData& trigData) const
2702 {
2703   // see header file for function documentation
2704   if (!fpCTPData) {
2705     static bool bWarningThrown=false;
2706     if (!bWarningThrown) HLTError("Trigger classes not initialized, use SetupCTPData from DoInit()");
2707     bWarningThrown=true;
2708     return false;
2709   }
2710
2711   return fpCTPData->EvaluateCTPTriggerClass(expression, trigData);
2712 }
2713
2714 int AliHLTComponent::CheckCTPTrigger(const char* name) const
2715 {
2716   // see header file for function documentation
2717   if (!fpCTPData) {
2718     static bool bWarningThrown=false;
2719     if (!bWarningThrown) HLTError("Trigger classes not initialized, use SetupCTPData from DoInit()");
2720     bWarningThrown=true;
2721     return false;
2722   }
2723
2724   return fpCTPData->CheckTrigger(name);
2725 }
2726
2727 Double_t AliHLTComponent::GetBz()
2728 {
2729   // Returns Bz.
2730   return AliHLTMisc::Instance().GetBz();
2731 }
2732
2733 Double_t AliHLTComponent::GetBz(const Double_t *r)
2734 {
2735   // Returns Bz (kG) at the point "r" .
2736   return AliHLTMisc::Instance().GetBz(r);
2737 }
2738
2739 void AliHLTComponent::GetBxByBz(const Double_t r[3], Double_t b[3])
2740 {
2741   // Returns Bx, By and Bz (kG) at the point "r" .
2742   AliHLTMisc::Instance().GetBxByBz(r, b);
2743 }