]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/sim/AliHLTOUTComponent.cxx
added new function to extract DDL id from HLT origin and specification, code cleanup
[u/mrichter/AliRoot.git] / HLT / sim / AliHLTOUTComponent.cxx
1 // $Id$
2
3 //**************************************************************************
4 //* This file is property of and copyright by the ALICE HLT Project        * 
5 //* ALICE Experiment at CERN, All rights reserved.                         *
6 //*                                                                        *
7 //* Primary Authors: Matthias Richter <Matthias.Richter@ift.uib.no>        *
8 //*                  for The ALICE HLT Project.                            *
9 //*                                                                        *
10 //* Permission to use, copy, modify and distribute this software and its   *
11 //* documentation strictly for non-commercial purposes is hereby granted   *
12 //* without fee, provided that the above copyright notice appears in all   *
13 //* copies and that both the copyright notice and this permission notice   *
14 //* appear in the supporting documentation. The authors make no claims     *
15 //* about the suitability of this software for any purpose. It is          *
16 //* provided "as is" without express or implied warranty.                  *
17 //**************************************************************************
18
19 //  @file   AliHLTOUTComponent.cxx
20 //  @author Matthias Richter
21 //  @date   
22 //  @brief  The HLTOUT data sink component similar to HLTOUT nodes
23 //  @note   Used in the AliRoot environment only.
24
25 #if __GNUC__>= 3
26 using namespace std;
27 #endif
28
29 #include <cassert>
30 //#include <iostream>
31 #include "AliHLTOUTComponent.h"
32 #include "AliHLTOUT.h"
33 #include "AliHLTHOMERLibManager.h"
34 #include "AliHLTHOMERWriter.h"
35 #include "AliHLTErrorGuard.h"
36 #include "AliDAQ.h" // equipment Ids
37 #include "AliRawDataHeader.h" // Common Data Header 
38 #include <TDatime.h> // seed for TRandom
39 #include <TRandom.h> // random int generation for DDL no
40 #include <TFile.h>
41 #include <TTree.h>
42 #include <TArrayC.h>
43 #include <TSystem.h>
44
45 /** ROOT macro for the implementation of ROOT specific class methods */
46 ClassImp(AliHLTOUTComponent)
47
48 AliHLTOUTComponent::AliHLTOUTComponent(EType type)
49   : AliHLTOfflineDataSink()
50   , fWriters()
51   , fNofDDLs(10)
52   , fIdFirstDDL(7680) // 0x1e<<8
53   , fBuffer()
54   , fpLibManager(NULL)
55   , fOptions(0)
56   , fDigitFileName("HLT.Digits.root")
57   , fpDigitFile(NULL)
58   , fpDigitTree(NULL)
59   , fppDigitArrays(NULL)
60   , fReservedWriter(-1)
61   , fReservedData(0)
62   , fType(type)
63 {
64   // see header file for class documentation
65   // or
66   // refer to README to build package
67   // or
68   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
69
70   // I guess DDL definitions should never change any more
71   assert(fNofDDLs==AliDAQ::NumberOfDdls("HLT"));
72   fNofDDLs=AliDAQ::NumberOfDdls("HLT");
73   
74   /* AliDAQ::DdlIDOffset returns wrong offset for HLT links
75   assert(fIdFirstDDL==AliDAQ::DdlIDOffset("HLT"));
76   fIdFirstDDL=AliDAQ::DdlIDOffset("HLT");
77   */
78
79   if (fType!=kGlobal && fType!=kDigits && fType!=kRaw) {
80     ALIHLTERRORGUARD(1, "invalid component type %d", fType);
81   }
82 }
83
84 int AliHLTOUTComponent::fgOptions=kWriteRawFiles|kWriteDigits;
85
86 AliHLTOUTComponent::~AliHLTOUTComponent()
87 {
88   // see header file for class documentation
89   if (fpLibManager) delete fpLibManager;
90   fpLibManager=NULL;
91 }
92
93 const char* AliHLTOUTComponent::GetComponentID()
94 {
95   // see header file for class documentation
96   switch (fType) {
97   case kDigits: return "HLTOUTdigits";
98   case kRaw:    return "HLTOUTraw";
99   case kGlobal:
100   default:
101     return "HLTOUT";
102   }
103   return NULL;
104 }
105
106 void AliHLTOUTComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list)
107 {
108   // see header file for class documentation
109   list.clear();
110   list.push_back(kAliHLTAnyDataType);
111 }
112
113 AliHLTComponent* AliHLTOUTComponent::Spawn()
114 {
115   // see header file for class documentation
116   return new AliHLTOUTComponent(fType);
117 }
118
119 int AliHLTOUTComponent::DoInit( int argc, const char** argv )
120 {
121   // see header file for class documentation
122   int iResult=0;
123
124   switch (fType) {
125   case kDigits:
126     fOptions|=kWriteDigits; fOptions&=~kWriteRawFiles;
127     HLTInfo("initializing HLTOUT component for digits generation");
128     break;
129   case kRaw:
130     fOptions|=kWriteRawFiles; fOptions&=~kWriteDigits;
131     HLTInfo("initializing HLTOUT component for raw data generation");
132     break;
133   case kGlobal:
134   default:
135     fOptions=fgOptions;
136   }
137
138   if ((iResult=ConfigureFromArgumentString(argc, argv))<0) return iResult;
139
140   // Make sure there is no library manager before we try and create a new one.
141   if (fpLibManager) {
142     delete fpLibManager;
143     fpLibManager=NULL;
144   }
145   
146   // Create a new library manager and allocate the appropriate number of
147   // HOMER writers for the HLTOUT component.
148   fpLibManager=new AliHLTHOMERLibManager;
149   if (fpLibManager) {
150     int writerNo=0;
151     for (writerNo=0; writerNo<fNofDDLs; writerNo++) {
152       AliHLTMonitoringWriter* pWriter=fpLibManager->OpenWriter();
153       if (pWriter) {
154         HLTDebug("HOMER writer %p added", pWriter);
155         fWriters.push_back(pWriter);
156       } else {
157         HLTError("can not open HOMER writer");
158         iResult=-ENODEV;
159         break;
160       }
161     }
162   } else {
163     iResult=-ENOMEM;
164   }
165
166   return iResult;
167 }
168
169 int AliHLTOUTComponent::ScanConfigurationArgument(int argc, const char** argv)
170 {
171   // see header file for class documentation
172   if (argc<=0) return 0;
173   int i=0;
174   TString argument=argv[i];
175   const char* key="";
176
177   // -links n
178   // specify number of ddl links
179   if (argument.CompareTo("-links")==0) {
180     if (++i>=argc) return -EPROTO;
181     TString parameter(argv[i]);
182     parameter.Remove(TString::kLeading, ' '); // remove all blanks
183     if (parameter.IsDigit()) {
184       fNofDDLs=parameter.Atoi();
185     } else {
186       HLTError("wrong parameter for argument %s, number expected", argument.Data());
187       return -EINVAL;
188     }
189
190     return 2;
191   } 
192
193   // -digitfile name
194   if (argument.CompareTo("-digitfile")==0) {
195     if (++i>=argc) return -EPROTO;
196     fDigitFileName=argv[i];
197
198     return 2;
199   }
200
201   // -rawout
202   key="-rawout";
203   if (argument.Contains(key)) {
204     argument.ReplaceAll(key, "");
205     if (argument.IsNull()) {
206       fOptions|=kWriteRawFiles;
207     } else if (argument.CompareTo("=off")==0) {
208       fOptions&=~kWriteRawFiles;
209     } else if (argument.CompareTo("=on")==0) {
210       fOptions|=kWriteRawFiles;
211     } else {
212       HLTError("invalid parameter for argument %s: possible %s=off/%s=on", key, key, key);
213       return -EPROTO;
214     }
215
216     return 1;
217   }
218
219   // -digitout
220   key="-digitout";
221   if (argument.Contains(key)) {
222     argument.ReplaceAll(key, "");
223     if (argument.IsNull()) {
224       fOptions|=kWriteDigits;
225     } else if (argument.CompareTo("=off")==0) {
226       fOptions&=~kWriteDigits;
227     } else if (argument.CompareTo("=on")==0) {
228       fOptions|=kWriteDigits;
229     } else {
230       HLTError("invalid parameter for argument %s: possible %s=off/%s=on", key, key, key);
231       return -EPROTO;
232     }
233
234     return 1;
235   }
236
237   // unknown argument
238   return -EINVAL;
239 }
240
241 int AliHLTOUTComponent::DoDeinit()
242 {
243   // see header file for class documentation
244   int iResult=0;
245
246   if (fpLibManager) {
247     AliHLTMonitoringWriterPVector::iterator element=fWriters.begin();
248     while (element!= fWriters.end()) {
249       assert(*element);
250       // wanted to have a dynamic_cast<AliHLTHOMERWriter*> here, but this results into
251       // undefined symbol when loading the library
252       if (*element!=NULL) {
253         (*element)->Clear();
254         fpLibManager->DeleteWriter((AliHLTHOMERWriter*)(*element));
255       } else {
256         HLTError("writer instance is NULL");
257       }
258       element=fWriters.erase(element);
259     }
260   }
261   if (fpLibManager) {
262     delete fpLibManager;
263     fpLibManager=NULL;
264   }
265
266   if (fpDigitTree) {
267     delete fpDigitTree;
268     fpDigitTree=NULL;
269   }
270
271   if (fpDigitFile) {
272     fpDigitFile->Close();
273     delete fpDigitFile;
274     fpDigitFile=NULL;
275   }
276
277   if (fppDigitArrays) {
278     for (int i=0; i<fNofDDLs; i++) {
279       if (fppDigitArrays[i]) delete fppDigitArrays[i];
280     }
281     delete[] fppDigitArrays;
282     fppDigitArrays=NULL;
283   }
284   
285   return iResult;
286 }
287
288 int AliHLTOUTComponent::DumpEvent( const AliHLTComponentEventData& evtData,
289                          const AliHLTComponentBlockData* blocks, 
290                          AliHLTComponentTriggerData& /*trigData*/ )
291 {
292   // see header file for class documentation
293   int iResult=0;
294   HLTInfo("write %d output block(s)", evtData.fBlockCnt);
295   int writerNo=0;
296   int blockCount=0;
297   AliHLTUInt32_t eventType=gkAliEventTypeUnknown;
298   bool bIsDataEvent=IsDataEvent(&eventType);
299   if (iResult>=0) {
300     homer_uint64 homerHeader[kCount_64b_Words];
301     HOMERBlockDescriptor homerDescriptor(homerHeader);
302     for (int n=0; n<(int)evtData.fBlockCnt; n++ ) {
303       if (blocks[n].fDataType==kAliHLTDataTypeEvent ||
304           blocks[n].fDataType==kAliHLTDataTypeSOR ||
305           blocks[n].fDataType==kAliHLTDataTypeEOR ||
306           blocks[n].fDataType==kAliHLTDataTypeComConf ||
307           blocks[n].fDataType==kAliHLTDataTypeUpdtDCS)
308         {
309           // the special events have to be ignored.
310           continue;
311         }
312       if (!bIsDataEvent &&
313           (blocks[n].fDataType!=kAliHLTDataTypeComponentTable))
314         {
315           // In simulation, there are no SOR and EOR events created. Thats
316           // why all data blocks of those events are currently ignored.
317           // Strictly speaking, components should not create output blocks
318           // on the SOR/EOR event
319           //
320           // Exeptions: some blocks are added, the buffer must be prepared and
321           // kept since the pointers will be invalid
322           // - kAliHLTDataTypeComponentTable component table entries
323           continue;
324         }
325       memset( homerHeader, 0, sizeof(homer_uint64)*kCount_64b_Words );
326       homerDescriptor.Initialize();
327       // for some traditional reason the TCPDumpSubscriber swaps the bytes
328       // of the data type id and data type origin. Actually I do not understand
329       // the corresponding code line
330       // homerBlock.SetType( blocks[n].fDataType.fID );
331       // this compiles in the PubSub framework and in addition does a byte swap
332       homer_uint64 id=0;
333       homer_uint64 origin=0;
334       memcpy(&id, blocks[n].fDataType.fID, sizeof(homer_uint64));
335       memcpy(((AliHLTUInt8_t*)&origin)+sizeof(homer_uint32), blocks[n].fDataType.fOrigin, sizeof(homer_uint32));
336       homerDescriptor.SetType(AliHLTOUT::ByteSwap64(id));
337       homerDescriptor.SetSubType1(AliHLTOUT::ByteSwap64(origin));
338       homerDescriptor.SetSubType2(blocks[n].fSpecification);
339       homerDescriptor.SetBlockSize(blocks[n].fSize);
340       if (bIsDataEvent) {
341         writerNo=ShuffleWriters(fWriters, blocks[n].fSize);
342       }
343       assert(writerNo>=0 && writerNo<(int)fWriters.size());
344       // I'm puzzled by the different headers, buffers etc. used in the
345       // HOMER writer/data. In additional, there is no type check as there
346       // are void pointers used and names mixed.
347       // It seems that HOMERBlockDescriptor is just a tool to set the
348       // different fields in the homer header, which is an array of 64 bit
349       // words.
350       fWriters[writerNo]->AddBlock(homerHeader, blocks[n].fPtr);
351       blockCount++;
352     }
353   }
354
355   if (iResult>=0 && !bIsDataEvent && fNofDDLs>=2) {
356     // data blocks from a special event are kept to be added to the
357     // following event. In the current implementation at least 2 DDLs
358     // are required to allow to keep the blocks of the SOR event and
359     // include it in the first event. If only one writer is available
360     // the blocks are ignored. For the moment this is not expexted to
361     // be a problem since components should not gererate anything on
362     // SOR/EOR. The only case is the list of AliHLTComponentTableEntry
363     // transmitted for component statistics in debug mode.
364     if (fReservedWriter>=0) {
365       HLTWarning("overriding previous buffer of non-data event data blocks");
366     }
367     const AliHLTUInt8_t* pBuffer=NULL;
368     int bufferSize=0;
369     // TODO: not yet clear whether it is smart to send the event id of
370     // this special event or if it should be set from the id of the
371     // following event where the data will be added
372     if (blockCount>0 && (bufferSize=FillOutputBuffer(evtData.fEventID, fWriters[writerNo], pBuffer))>0) {
373       fReservedWriter=writerNo;
374       fReservedData=bufferSize;
375     }
376     fWriters[writerNo]->Clear();
377   } else if (iResult>=0 && !bIsDataEvent && fNofDDLs<2 && blockCount>0) {
378     HLTWarning("ignoring %d block(s) for special event of type %d: at least 2 DDLs are required", blockCount, eventType);
379   }
380
381   if (iResult>=0 && bIsDataEvent) {
382     iResult=Write(GetEventCount(), GetRunLoader());
383   }
384
385   return iResult;
386 }
387
388
389 int AliHLTOUTComponent::FillESD(int /*eventNo*/, AliRunLoader* /*runLoader*/, AliESDEvent* /*esd*/)
390 {
391   // see header file for class documentation
392   // 2010-04-14 nothing to do any more. The data is written at the end of
393   // DumpEvent
394   return 0;
395 }
396
397 int AliHLTOUTComponent::Write(int eventNo, AliRunLoader* runLoader)
398 {
399   // see header file for class documentation
400   int iResult=0;
401
402   if (fWriters.size()==0) return 0;
403
404   if (fReservedWriter>=0) {
405     if (fOptions&kWriteDigits) WriteDigitArray(fReservedWriter, &fBuffer[0], fReservedData);
406     if (fOptions&kWriteRawFiles) WriteRawFile(eventNo, runLoader, fReservedWriter, &fBuffer[0], fReservedData);
407     fReservedData=0;
408   }
409
410   // search for the writer with the biggest data volume in order to allocate the
411   // output buffer of sufficient size
412   vector<int> sorted;
413   for (size_t i=0; i<fWriters.size(); i++) {
414     if ((int)i==fReservedWriter) continue;    
415     assert(fWriters[i]);
416     if (fWriters[i]) {
417       if (sorted.size()==0 || fWriters[i]->GetTotalMemorySize()<=fWriters[sorted[0]]->GetTotalMemorySize()) {
418         sorted.push_back(i);
419       } else {
420         sorted.insert(sorted.begin(), i);
421       }
422     }
423   }
424   fReservedWriter=-1;
425
426   vector<int>::iterator ddlno=sorted.begin();
427   while (ddlno!=sorted.end()) {
428     const AliHLTUInt8_t* pBuffer=NULL;
429     int bufferSize=0;
430     
431     if ((bufferSize=FillOutputBuffer(eventNo, fWriters[*ddlno], pBuffer))>0) {
432       if (fOptions&kWriteDigits) WriteDigitArray(*ddlno, pBuffer, bufferSize);
433       if (fOptions&kWriteRawFiles) WriteRawFile(eventNo, runLoader, *ddlno, pBuffer, bufferSize);
434     }
435     fWriters[*ddlno]->Clear();
436     ddlno++;
437   }
438   if (fOptions&kWriteDigits) WriteDigits(eventNo, runLoader);
439   return iResult;
440 }
441
442 int AliHLTOUTComponent::ShuffleWriters(AliHLTMonitoringWriterPVector &list, AliHLTUInt32_t /*size*/)
443 {
444   // see header file for class documentation
445   int iResult=-ENOENT;
446   assert(list.size()>0);
447   if (list.size()==0) return iResult;
448   vector<int> writers;
449   size_t i=0;
450   for (i=0; i<list.size(); i++) {
451     if ((int)i==fReservedWriter) continue;
452     if (list[i]->GetTotalMemorySize()==0)
453       writers.push_back(i);
454     else if (iResult<0 ||
455              list[i]->GetTotalMemorySize()<list[iResult]->GetTotalMemorySize())
456       iResult=i;
457       
458   }
459   if (writers.size()>0) {
460     iResult=writers[0];
461     if (writers.size()>0) {
462       // shuffle among the empty writers
463       TDatime dt;
464       TRandom rand;
465       rand.SetSeed(dt.Get()*(iResult+1));
466       i=rand.Integer(writers.size()-1);
467       assert(i>0 && i<writers.size()-1);
468       iResult=writers[i];
469     }
470   } else {
471     // take the writer with the least data volume
472     assert(iResult>=0);
473   }
474   return iResult;
475 }
476
477 int AliHLTOUTComponent::FillOutputBuffer(int eventNo, AliHLTMonitoringWriter* pWriter, const AliHLTUInt8_t* &pBuffer)
478 {
479   // see header file for class documentation
480   int iResult=0;
481   unsigned int bufferSize=0;
482
483   // space for common data header
484   bufferSize+=sizeof(AliRawDataHeader);
485   assert(sizeof(AliRawDataHeader)==32);
486
487   // space for HLT event header
488   bufferSize+=sizeof(AliHLTOUT::AliHLTOUTEventHeader);
489
490   // space for payload from the writer
491   if (pWriter) bufferSize+=pWriter->GetTotalMemorySize();
492
493   // payload data must be aligned to 32bit
494   bufferSize=(bufferSize+3)/4;
495   bufferSize*=4;
496
497   if (bufferSize>fBuffer.size())
498     fBuffer.resize(bufferSize);
499
500   // reset the last 32bit word, rest will be overwritten
501   memset(&fBuffer[bufferSize-4], 0, 4);
502
503   if (bufferSize<=fBuffer.size()) {
504     AliRawDataHeader* pCDH=reinterpret_cast<AliRawDataHeader*>(&fBuffer[0]);
505     AliHLTOUT::AliHLTOUTEventHeader* pHLTH=reinterpret_cast<AliHLTOUT::AliHLTOUTEventHeader*>(&fBuffer[sizeof(AliRawDataHeader)]);
506     *pCDH = AliRawDataHeader();  // Fill with default values.
507     memset(pHLTH, 0, sizeof(AliHLTOUT::AliHLTOUTEventHeader));
508
509     if (pWriter) {
510       // copy payload
511       pWriter->Copy(&fBuffer[sizeof(AliRawDataHeader)+sizeof(AliHLTOUT::AliHLTOUTEventHeader)], 0, 0, 0, 0);
512       pHLTH->fLength=pWriter->GetTotalMemorySize();
513       // set status bit to indicate HLT payload
514       pCDH->fStatusMiniEventID|=0x1<<(AliHLTOUT::kCDHStatusFlagsOffset+AliHLTOUT::kCDHFlagsHLTPayload);
515     }
516     pHLTH->fLength+=sizeof(AliHLTOUT::AliHLTOUTEventHeader);
517     // pHLTH->fEventIDLow is already set to zero in memset above.
518     pHLTH->fEventIDLow = eventNo;
519     // version does not really matter since we do not add decision data
520     pHLTH->fVersion=AliHLTOUT::kVersion1;
521
522     pCDH->fSize=bufferSize;
523     pCDH->fStatusMiniEventID|=0x1<<(AliHLTOUT::kCDHStatusFlagsOffset + AliHLTOUT::kCDHFlagsHLTPayload);
524     
525     pBuffer=&fBuffer[0];
526     iResult=(int)bufferSize;
527   } else {
528     pBuffer=NULL;
529     iResult=-ENOMEM;
530   }
531
532   return iResult;
533 }
534
535 int AliHLTOUTComponent::WriteDigitArray(int hltddl, const AliHLTUInt8_t* pBuffer, unsigned int bufferSize)
536 {
537   // see header file for class documentation
538   int iResult=0;
539   assert(hltddl<fNofDDLs);
540   if (hltddl>=fNofDDLs) return -ERANGE;
541
542   if (!fppDigitArrays) {
543     fppDigitArrays=new TArrayC*[fNofDDLs];
544     if (fppDigitArrays) {
545       for (int i=0; i<fNofDDLs; i++) {
546         fppDigitArrays[i]=new TArrayC(0);
547       }
548     }
549   }
550   if (fppDigitArrays && fppDigitArrays[hltddl]) {
551     fppDigitArrays[hltddl]->Set(bufferSize, reinterpret_cast<const Char_t*>(pBuffer));
552   } else {
553     iResult=-ENOMEM;    
554   }
555   return iResult;
556 }
557
558 int AliHLTOUTComponent::WriteDigits(int /*eventNo*/, AliRunLoader* /*runLoader*/)
559 {
560   // see header file for class documentation
561   int iResult=0;
562   if (!fpDigitFile) {
563     fpDigitFile=new TFile(fDigitFileName, "RECREATE");
564   }
565   if (fpDigitFile && !fpDigitFile->IsZombie()) {
566     if (!fpDigitTree) {
567       fpDigitTree=new TTree("rawhltout","HLTOUT raw data");
568       if (fpDigitTree && fppDigitArrays) {
569         for (int i=0; i<fNofDDLs; i++) {
570           const char* branchName=AliDAQ::DdlFileName("HLT", i);
571           if (fppDigitArrays[i]) fpDigitTree->Branch(branchName, "TArrayC", &fppDigitArrays[i], 32000/*just as the default*/, 0);
572         }
573       }
574     }
575     if (fpDigitTree) {
576 #ifdef __DEBUG
577       int res=fpDigitTree->Fill();
578       HLTDebug("writing digit tree: %d", res);
579       fpDigitFile->cd();
580       res=fpDigitTree->Write("",TObject::kOverwrite);
581       HLTDebug("writing digit tree: %d", res);
582 #else
583       fpDigitTree->Fill();
584       fpDigitFile->cd();
585       fpDigitTree->Write("",TObject::kOverwrite);
586 #endif
587       if (fppDigitArrays) for (int i=0; i<fNofDDLs; i++) {
588         if (fppDigitArrays[i]) fppDigitArrays[i]->Set(0);
589       }
590     }
591   } else {
592     const char* errorMsg="";
593     if (GetEventCount()==5) {
594       errorMsg=" (suppressing further error messages)";
595     }
596     if (GetEventCount()<5) {
597       HLTError("can not open HLT digit file %s%s", fDigitFileName.Data(), errorMsg);
598     }
599     iResult=-EBADF;
600   }
601   return iResult;
602 }
603
604 int AliHLTOUTComponent::WriteRawFile(int eventNo, AliRunLoader* /*runLoader*/, int hltddl, const AliHLTUInt8_t* pBuffer, unsigned int bufferSize)
605 {
606   // see header file for class documentation
607   int iResult=0;
608   const char* fileName=AliDAQ::DdlFileName("HLT", hltddl);
609   assert(fileName!=NULL);
610   TString filePath;
611   filePath.Form("raw%d/", eventNo);
612   if (gSystem->AccessPathName(filePath)) {
613     TString command="mkdir "; command+=filePath;
614     gSystem->Exec(command);
615   }
616   filePath+=fileName;
617   if (fileName) {
618     ios::openmode filemode=(ios::openmode)0;
619     ofstream rawfile(filePath.Data(), filemode);
620     if (rawfile.good()) {
621       if (pBuffer && bufferSize>0) {
622         rawfile.write(reinterpret_cast<const char*>(pBuffer), bufferSize);
623       } else {
624         HLTWarning("writing zero length raw data file %s");
625       }
626       HLTDebug("wrote %d byte(s) to file %s", bufferSize, filePath.Data());
627     } else {
628       HLTError("can not open file %s for writing", filePath.Data());
629       iResult=-EBADF;
630     }
631     rawfile.close();
632   }
633   return iResult;
634 }
635
636 void AliHLTOUTComponent::SetGlobalOption(unsigned int options)
637 {
638   // see header file for class documentation
639   fgOptions|=options;
640 }
641
642 void AliHLTOUTComponent::ClearGlobalOption(unsigned int options)
643 {
644   // see header file for class documentation
645   fgOptions&=~options;
646 }
647
648 bool AliHLTOUTComponent::TestGlobalOption(unsigned int option)
649 {
650   // see header file for class documentation
651   return (fgOptions&option)!=0;
652 }