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