]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/sim/AliHLTOUTComponent.cxx
f4ebcadbee862994ac2d20c46b02fe14b834dc1a
[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 */
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 "AliDAQ.h" // equipment Ids
36 #include "AliRawDataHeader.h" // Common Data Header 
37 #include <TDatime.h> // seed for TRandom
38 #include <TRandom.h> // random int generation for DDL no
39 #include <TFile.h>
40 #include <TTree.h>
41 #include <TArrayC.h>
42
43 /** ROOT macro for the implementation of ROOT specific class methods */
44 ClassImp(AliHLTOUTComponent)
45
46 AliHLTOUTComponent::AliHLTOUTComponent()
47   :
48   AliHLTOfflineDataSink(),
49   fWriters(),
50   fNofDDLs(10),
51   fIdFirstDDL(7680), // 0x1e<<8
52   fBuffer(),
53   fpLibManager(NULL),
54   fpDigitFile(NULL),
55   fpDigitTree(NULL),
56   fppDigitArrays(NULL)
57 {
58   // see header file for class documentation
59   // or
60   // refer to README to build package
61   // or
62   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
63
64   // I guess DDL definitions should never change any more
65   assert(fNofDDLs==AliDAQ::NumberOfDdls("HLT"));
66   fNofDDLs=AliDAQ::NumberOfDdls("HLT");
67   
68   /* AliDAQ::DdlIDOffset returns wrong offset for HLT links
69   assert(fIdFirstDDL==AliDAQ::DdlIDOffset("HLT"));
70   fIdFirstDDL=AliDAQ::DdlIDOffset("HLT");
71   */
72 }
73
74 int AliHLTOUTComponent::fgOptions=kWriteRawFiles|kWriteDigits;
75
76 AliHLTOUTComponent::~AliHLTOUTComponent()
77 {
78   // see header file for class documentation
79   if (fpLibManager) delete fpLibManager;
80   fpLibManager=NULL;
81 }
82
83 const char* AliHLTOUTComponent::GetComponentID()
84 {
85   // see header file for class documentation
86   return "HLTOUT";
87 }
88
89 void AliHLTOUTComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list)
90 {
91   // see header file for class documentation
92   list.clear();
93   list.push_back(kAliHLTAnyDataType);
94 }
95
96 AliHLTComponent* AliHLTOUTComponent::Spawn()
97 {
98   // see header file for class documentation
99   return new AliHLTOUTComponent;
100 }
101
102 int AliHLTOUTComponent::DoInit( int argc, const char** argv )
103 {
104   // see header file for class documentation
105   int iResult=0;
106   TString argument="";
107   int bMissingParam=0;
108   for (int i=0; i<argc && iResult>=0; i++) {
109     argument=argv[i];
110     if (argument.IsNull()) continue;
111
112     // -links
113     if (argument.CompareTo("-links")==0) {
114       if ((bMissingParam=(++i>=argc))) break;
115       TString parameter(argv[i]);
116       parameter.Remove(TString::kLeading, ' '); // remove all blanks
117       if (parameter.IsDigit()) {
118         fNofDDLs=parameter.Atoi();
119       } else {
120         HLTError("wrong parameter for argument %s, number expected", argument.Data());
121         iResult=-EINVAL;
122       }
123     } else {
124       HLTError("unknown argument %s", argument.Data());
125       iResult=-EINVAL;
126       break;
127     }
128   }
129   if (bMissingParam) {
130     HLTError("missing parameter for argument %s", argument.Data());
131     iResult=-EINVAL;
132   }
133   if (iResult>=0) {
134   }
135
136   fpLibManager=new AliHLTHOMERLibManager;
137   if (fpLibManager) {
138     int writerNo=0;
139     for (writerNo=0; writerNo<fNofDDLs; writerNo++) {
140       AliHLTMonitoringWriter* pWriter=fpLibManager->OpenWriter();
141       if (pWriter) {
142         HLTDebug("HOMER writer %p added", pWriter);
143         fWriters.push_back(pWriter);
144       } else {
145         HLTError("can not open HOMER writer");
146         iResult=-ENODEV;
147         break;
148       }
149     }
150   } else {
151     iResult=-ENOMEM;
152   }
153
154   return iResult;
155 }
156
157 int AliHLTOUTComponent::DoDeinit()
158 {
159   // see header file for class documentation
160   int iResult=0;
161
162   if (fpLibManager) {
163     AliHLTMonitoringWriterPVector::iterator element=fWriters.begin();
164     while (element!= fWriters.end()) {
165       assert(*element);
166       // wanted to have a dynamic_cast<AliHLTHOMERWriter*> here, but this results into
167       // undefined symbol when loading the library
168       (*element)->Clear();
169       if (*element!=NULL) fpLibManager->DeleteWriter((AliHLTHOMERWriter*)(*element));
170       element=fWriters.erase(element);
171     }
172   }
173
174   if (fpDigitTree) {
175     delete fpDigitTree;
176     fpDigitTree=NULL;
177   }
178
179   if (fpDigitFile) {
180     fpDigitFile->Close();
181     delete fpDigitFile;
182     fpDigitFile=NULL;
183   }
184
185   if (fppDigitArrays) {
186     for (int i=0; i<fNofDDLs; i++) {
187       if (fppDigitArrays[i]) delete fppDigitArrays[i];
188     }
189     delete[] fppDigitArrays;
190     fppDigitArrays=NULL;
191   }
192   
193   return iResult;
194 }
195
196 int AliHLTOUTComponent::DumpEvent( const AliHLTComponentEventData& evtData,
197                          const AliHLTComponentBlockData* blocks, 
198                          AliHLTComponentTriggerData& /*trigData*/ )
199 {
200   // see header file for class documentation
201   int iResult=0;
202   HLTInfo("write %d output block(s)", evtData.fBlockCnt);
203   if (iResult>=0) {
204     homer_uint64 homerHeader[kCount_64b_Words];
205     HOMERBlockDescriptor homerDescriptor(homerHeader);
206     for (int n=0; n<(int)evtData.fBlockCnt; n++ ) {
207       memset( homerHeader, 0, sizeof(homer_uint64)*kCount_64b_Words );
208       homerDescriptor.Initialize();
209       // for some traditional reason the TCPDumpSubscriber swaps the bytes
210       // of the data type id and data type origin. Actually I do not understand
211       // the corresponding code line
212       // homerBlock.SetType( blocks[n].fDataType.fID );
213       // this compiles in the PubSub framework and in addition does a byte swap
214       homer_uint64 id=0;
215       homer_uint64 origin=0;
216       memcpy(&id, blocks[n].fDataType.fID, sizeof(homer_uint64));
217       memcpy(((AliHLTUInt8_t*)&origin)+sizeof(homer_uint32), blocks[n].fDataType.fOrigin, sizeof(homer_uint32));
218       homerDescriptor.SetType(AliHLTOUT::ByteSwap64(id));
219       homerDescriptor.SetSubType1(AliHLTOUT::ByteSwap64(origin));
220       homerDescriptor.SetSubType2(blocks[n].fSpecification);
221       homerDescriptor.SetBlockSize(blocks[n].fSize);
222       int writerNo=ShuffleWriters(fWriters, blocks[n].fSize);
223       assert(writerNo>=0 && writerNo<(int)fWriters.size());
224       // I'm puzzled by the different headers, buffers etc. used in the
225       // HOMER writer/data. In additional, there is no type check as there
226       // are void pointers used and names mixed.
227       // It seems that HOMERBlockDescriptor is just a tool to set the
228       // different fields in the homer header, which is an array of 64 bit
229       // words.
230       fWriters[writerNo]->AddBlock(homerHeader, blocks[n].fPtr);
231     }
232   }
233
234   return iResult;
235 }
236
237 int AliHLTOUTComponent::FillESD(int eventNo, AliRunLoader* runLoader, AliESDEvent* /*esd*/)
238 {
239   // see header file for class documentation
240   int iResult=0;
241   if (fWriters.size()==0) return 0;
242   
243   // search for the writer with the biggest data volume in order to allocate the
244   // output buffer of sufficient size
245   vector<int> sorted;
246   for (size_t i=0; i<fWriters.size(); i++) {
247     assert(fWriters[i]);
248     if (fWriters[i]) {
249       if (sorted.size()==0 || fWriters[i]->GetTotalMemorySize()<=fWriters[sorted[0]]->GetTotalMemorySize()) {
250         sorted.push_back(i);
251       } else {
252         sorted.insert(sorted.begin(), i);
253       }
254     }
255   }
256
257   vector<int>::iterator ddlno=sorted.begin();
258   while (ddlno!=sorted.end()) {
259     const AliHLTUInt8_t* pBuffer=NULL;
260     int bufferSize=0;
261     
262     if ((bufferSize=FillOutputBuffer(eventNo, fWriters[*ddlno], pBuffer))>0) {
263       if (fgOptions&kWriteDigits) WriteDigitArray(*ddlno, pBuffer, bufferSize);
264       if (fgOptions&kWriteRawFiles) WriteRawFile(eventNo, runLoader, *ddlno, pBuffer, bufferSize);
265     }
266     fWriters[*ddlno]->Clear();
267     ddlno++;
268   }
269   if (fgOptions&kWriteDigits) WriteDigits(eventNo, runLoader);
270   return iResult;
271 }
272
273 int AliHLTOUTComponent::ShuffleWriters(AliHLTMonitoringWriterPVector &list, AliHLTUInt32_t /*size*/)
274 {
275   // see header file for class documentation
276   int iResult=-ENOENT;
277   assert(list.size()>0);
278   if (list.size()==0) return iResult;
279   vector<int> writers;
280   size_t i=0;
281   for (i=0; i<list.size(); i++) {
282     if (list[i]->GetTotalMemorySize()==0)
283       writers.push_back(i);
284     else if (iResult<0 ||
285              list[i]->GetTotalMemorySize()<list[iResult]->GetTotalMemorySize())
286       iResult=i;
287       
288   }
289   if (writers.size()>0) {
290     iResult=writers[0];
291     if (writers.size()>0) {
292       // shuffle among the empty writers
293       TDatime dt;
294       TRandom rand;
295       rand.SetSeed(dt.Get()*(iResult+1));
296       i=rand.Integer(writers.size()-1);
297       assert(i>0 && i<writers.size()-1);
298       iResult=writers[i];
299     }
300   } else {
301     // take the writer with the least data volume
302     assert(iResult>=0);
303   }
304   return iResult;
305 }
306
307 int AliHLTOUTComponent::FillOutputBuffer(int eventNo, AliHLTMonitoringWriter* pWriter, const AliHLTUInt8_t* &pBuffer)
308 {
309   // see header file for class documentation
310   int iResult=0;
311   unsigned int bufferSize=0;
312
313   // space for common data header
314   bufferSize+=sizeof(AliRawDataHeader);
315   assert(sizeof(AliRawDataHeader)==32);
316
317   // space for HLT event header
318   bufferSize+=sizeof(AliHLTOUT::AliHLTOUTEventHeader);
319
320   // space for payload from the writer
321   if (pWriter) bufferSize+=pWriter->GetTotalMemorySize();
322
323   if (bufferSize>fBuffer.size())
324     fBuffer.resize(bufferSize);
325
326   if (bufferSize<=fBuffer.size()) {
327     AliRawDataHeader* pCDH=reinterpret_cast<AliRawDataHeader*>(&fBuffer[0]);
328     AliHLTOUT::AliHLTOUTEventHeader* pHLTH=reinterpret_cast<AliHLTOUT::AliHLTOUTEventHeader*>(&fBuffer[sizeof(AliRawDataHeader)]);
329     *pCDH = AliRawDataHeader();  // Fill with default values.
330     memset(pHLTH, 0, sizeof(AliHLTOUT::AliHLTOUTEventHeader));
331
332     if (pWriter) {
333       // copy payload
334       pWriter->Copy(&fBuffer[sizeof(AliRawDataHeader)+sizeof(AliHLTOUT::AliHLTOUTEventHeader)], 0, 0, 0, 0);
335       pHLTH->fLength=pWriter->GetTotalMemorySize();
336       // set status bit to indicate HLT payload
337       pCDH->fStatusMiniEventID|=0x1<<(AliHLTOUT::kCDHStatusFlagsOffset+AliHLTOUT::kCDHFlagsHLTPayload);
338     }
339     pHLTH->fLength+=sizeof(AliHLTOUT::AliHLTOUTEventHeader);
340     // pHLTH->fEventIDLow is already set to zero in memset above.
341     pHLTH->fEventIDLow = eventNo;
342     // version does not really matter since we do not add decision data
343     pHLTH->fVersion=AliHLTOUT::kVersion1;
344
345     pCDH->fSize=sizeof(AliRawDataHeader)+pHLTH->fLength;
346     pCDH->fStatusMiniEventID|=0x1<<(AliHLTOUT::kCDHStatusFlagsOffset + AliHLTOUT::kCDHFlagsHLTPayload);
347     
348     pBuffer=&fBuffer[0];
349     iResult=(int)bufferSize;
350   } else {
351     pBuffer=NULL;
352     iResult=-ENOMEM;
353   }
354
355   return iResult;
356 }
357
358 int AliHLTOUTComponent::WriteDigitArray(int hltddl, const AliHLTUInt8_t* pBuffer, unsigned int bufferSize)
359 {
360   // see header file for class documentation
361   int iResult=0;
362   assert(hltddl<fNofDDLs);
363   if (hltddl>=fNofDDLs) return -ERANGE;
364
365   if (!fppDigitArrays) {
366     fppDigitArrays=new TArrayC*[fNofDDLs];
367     if (fppDigitArrays) {
368       for (int i=0; i<fNofDDLs; i++) {
369         fppDigitArrays[i]=new TArrayC(0);
370       }
371     }
372   }
373   if (fppDigitArrays && fppDigitArrays[hltddl]) {
374     fppDigitArrays[hltddl]->Set(bufferSize, reinterpret_cast<const Char_t*>(pBuffer));
375   } else {
376     iResult=-ENOMEM;    
377   }
378   return iResult;
379 }
380
381 int AliHLTOUTComponent::WriteDigits(int /*eventNo*/, AliRunLoader* /*runLoader*/)
382 {
383   // see header file for class documentation
384   int iResult=0;
385   const char* digitFileName="HLT.Digits.root";
386   if (!fpDigitFile) {
387     fpDigitFile=new TFile(digitFileName, "RECREATE");
388   }
389   if (fpDigitFile && !fpDigitFile->IsZombie()) {
390     if (!fpDigitTree) {
391       fpDigitTree=new TTree("rawhltout","HLTOUT raw data");
392       if (fpDigitTree && fppDigitArrays) {
393         for (int i=0; i<fNofDDLs; i++) {
394           const char* branchName=AliDAQ::DdlFileName("HLT", i);
395           if (fppDigitArrays[i]) fpDigitTree->Branch(branchName, "TArrayC", &fppDigitArrays[i], 32000/*just as the default*/, 0);
396         }
397       }
398     }
399     if (fpDigitTree) {
400       int res=fpDigitTree->Fill();
401       HLTDebug("writing digit tree: %d", res);
402       fpDigitFile->cd();
403       res=fpDigitTree->Write("",TObject::kOverwrite);
404       HLTDebug("writing digit tree: %d", res);
405       if (fppDigitArrays) for (int i=0; i<fNofDDLs; i++) {
406         if (fppDigitArrays[i]) fppDigitArrays[i]->Set(0);
407       }
408     }
409   } else {
410     const char* errorMsg="";
411     if (GetEventCount()==5) {
412       errorMsg=" (suppressing further error messages)";
413     }
414     if (GetEventCount()<5) {
415       HLTError("can not open HLT digit file %s%s", digitFileName, errorMsg);
416     }
417     iResult=-EBADF;
418   }
419   return iResult;
420 }
421
422 int AliHLTOUTComponent::WriteRawFile(int eventNo, AliRunLoader* /*runLoader*/, int hltddl, const AliHLTUInt8_t* pBuffer, unsigned int bufferSize)
423 {
424   // see header file for class documentation
425   int iResult=0;
426   const char* fileName=AliDAQ::DdlFileName("HLT", hltddl);
427   assert(fileName!=NULL);
428   TString filePath;
429   filePath.Form("raw%d/", eventNo);
430   filePath+=fileName;
431   if (fileName) {
432     ios::openmode filemode=(ios::openmode)0;
433     ofstream rawfile(filePath.Data(), filemode);
434     if (rawfile.good()) {
435       if (pBuffer && bufferSize>0) {
436         rawfile.write(reinterpret_cast<const char*>(pBuffer), bufferSize);
437       } else {
438         HLTWarning("writing zero length raw data file %s");
439       }
440       HLTDebug("wrote %d byte(s) to file %s", bufferSize, filePath.Data());
441     } else {
442       HLTError("can not open file %s for writing", filePath.Data());
443       iResult=-EBADF;
444     }
445     rawfile.close();
446   }
447   return iResult;
448 }
449
450 void AliHLTOUTComponent::SetGlobalOption(unsigned int options)
451 {
452   // see header file for class documentation
453   fgOptions|=options;
454 }
455
456 void AliHLTOUTComponent::ClearGlobalOption(unsigned int options)
457 {
458   // see header file for class documentation
459   fgOptions&=~options;
460 }