]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/sim/AliHLTOUTComponent.cxx
- handling of global AliHLTSystem singleton moved to BASE
[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 // see header file for class documentation
26 // or
27 // refer to README to build package
28 // or
29 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
30
31 #if __GNUC__>= 3
32 using namespace std;
33 #endif
34
35 #include <cassert>
36 //#include <iostream>
37 #include "AliHLTOUTComponent.h"
38 #include "AliHLTOUT.h"
39 #include "AliHLTHOMERLibManager.h"
40 #include "AliHLTHOMERWriter.h"
41 #include "AliDAQ.h" // equipment Ids
42 #include "AliRawDataHeader.h" // Common Data Header 
43 #include <TDatime.h> // seed for TRandom
44 #include <TRandom.h> // random int generation for DDL no
45
46 /** ROOT macro for the implementation of ROOT specific class methods */
47 ClassImp(AliHLTOUTComponent)
48
49 AliHLTOUTComponent::AliHLTOUTComponent()
50   :
51   AliHLTOfflineDataSink(),
52   fWriters(),
53   fNofDDLs(10),
54   fIdFirstDDL(7680), // 0x1e<<8
55   fBuffer(),
56   fpLibManager(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 nor 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   return iResult;
175 }
176
177 int AliHLTOUTComponent::DumpEvent( const AliHLTComponentEventData& evtData,
178                          const AliHLTComponentBlockData* blocks, 
179                          AliHLTComponentTriggerData& /*trigData*/ )
180 {
181   // see header file for class documentation
182   int iResult=0;
183   HLTInfo("write %d output block(s)", evtData.fBlockCnt);
184   if (iResult>=0) {
185     homer_uint64 homerHeader[kCount_64b_Words];
186     HOMERBlockDescriptor homerDescriptor(homerHeader);
187     for (int n=0; n<(int)evtData.fBlockCnt; n++ ) {
188       memset( homerHeader, 0, sizeof(homer_uint64)*kCount_64b_Words );
189       homerDescriptor.Initialize();
190       // for some traditional reason the TCPDumpSubscriber swaps the bytes
191       // of the data type id and data type origin. Actually I do not understand
192       // the corresponding code line
193       // homerBlock.SetType( blocks[n].fDataType.fID );
194       // this compiles in the PubSub framework and in addition does a byte swap
195       homer_uint64 id=0;
196       homer_uint64 origin=0;
197       memcpy(&id, blocks[n].fDataType.fID, sizeof(homer_uint64));
198       memcpy(((AliHLTUInt8_t*)&origin)+sizeof(homer_uint32), blocks[n].fDataType.fOrigin, sizeof(homer_uint32));
199       homerDescriptor.SetType(AliHLTOUT::ByteSwap64(id));
200       homerDescriptor.SetSubType1(AliHLTOUT::ByteSwap32(origin));
201       homerDescriptor.SetSubType2(blocks[n].fSpecification);
202       homerDescriptor.SetBlockSize(blocks[n].fSize);
203       int writerNo=ShuffleWriters(fWriters, blocks[n].fSize);
204       assert(writerNo>=0 && writerNo<(int)fWriters.size());
205       // I'm puzzled by the different headers, buffers etc. used in the
206       // HOMER writer/data. In additional, there is no type check as there
207       // are void pointers used and names mixed.
208       // It seems that HOMERBlockDescriptor is just a tool to set the
209       // different fields in the homer header, which is an array of 64 bit
210       // words.
211       fWriters[writerNo]->AddBlock(homerHeader, blocks[n].fPtr);
212     }
213   }
214
215   return iResult;
216 }
217
218 int AliHLTOUTComponent::FillESD(int eventNo, AliRunLoader* runLoader, AliESDEvent* /*esd*/)
219 {
220   // see header file for class documentation
221   int iResult=0;
222   if (fWriters.size()==0) return 0;
223   
224   // search for the writer with the biggest data volume in order to allocate the
225   // output buffer of sufficient size
226   vector<int> sorted;
227   for (size_t i=0; i<fWriters.size(); i++) {
228     assert(fWriters[i]);
229     if (fWriters[i]) {
230       if (sorted.size()==0 || fWriters[i]->GetTotalMemorySize()<=fWriters[sorted[0]]->GetTotalMemorySize()) {
231         sorted.push_back(i);
232       } else {
233         sorted.insert(sorted.begin(), i);
234       }
235     }
236   }
237
238   vector<int>::iterator ddlno=sorted.begin();
239   while (ddlno!=sorted.end()) {
240     const AliHLTUInt8_t* pBuffer=NULL;
241     int bufferSize=0;
242     
243     if ((bufferSize=FillOutputBuffer(eventNo, fWriters[*ddlno], pBuffer))>0) {
244       if (fgOptions&kWriteDigits) WriteDigits(eventNo, runLoader, *ddlno, pBuffer, bufferSize);
245       if (fgOptions&kWriteRawFiles) WriteRawFile(eventNo, runLoader, *ddlno, pBuffer, bufferSize);
246     }
247     fWriters[*ddlno]->Clear();
248     ddlno++;
249   }
250   return iResult;
251 }
252
253 int AliHLTOUTComponent::ShuffleWriters(AliHLTMonitoringWriterPVector &list, AliHLTUInt32_t /*size*/)
254 {
255   // see header file for class documentation
256   int iResult=-ENOENT;
257   assert(list.size()>0);
258   if (list.size()==0) return iResult;
259   vector<int> writers;
260   size_t i=0;
261   for (i=0; i<list.size(); i++) {
262     if (list[i]->GetTotalMemorySize()==0)
263       writers.push_back(i);
264     else if (iResult<0 ||
265              list[i]->GetTotalMemorySize()<list[iResult]->GetTotalMemorySize())
266       iResult=i;
267       
268   }
269   if (writers.size()>0) {
270     iResult=writers[0];
271     if (writers.size()>0) {
272       // shuffle among the empty writers
273       TDatime dt;
274       TRandom rand;
275       rand.SetSeed(dt.Get()*(iResult+1));
276       i=rand.Integer(writers.size()-1);
277       assert(i>0 && i<writers.size()-1);
278       iResult=writers[i];
279     }
280   } else {
281     // take the writer with the least data volume
282     assert(iResult>=0);
283   }
284   return iResult;
285 }
286
287 int AliHLTOUTComponent::FillOutputBuffer(int eventNo, AliHLTMonitoringWriter* pWriter, const AliHLTUInt8_t* &pBuffer)
288 {
289   // see header file for class documentation
290   int iResult=0;
291   unsigned int bufferSize=0;
292
293   // space for common data header
294   bufferSize+=sizeof(AliRawDataHeader);
295   assert(sizeof(AliRawDataHeader)==32);
296
297   // space for HLT event header
298   bufferSize+=sizeof(AliHLTOUT::AliHLTOUTEventHeader);
299
300   // space for payload from the writer
301   if (pWriter) bufferSize+=pWriter->GetTotalMemorySize();
302
303   if (bufferSize>fBuffer.size())
304     fBuffer.resize(bufferSize);
305
306   if (bufferSize<=fBuffer.size()) {
307     AliRawDataHeader* pCDH=reinterpret_cast<AliRawDataHeader*>(&fBuffer[0]);
308     AliHLTOUT::AliHLTOUTEventHeader* pHLTH=reinterpret_cast<AliHLTOUT::AliHLTOUTEventHeader*>(&fBuffer[sizeof(AliRawDataHeader)]);
309     memset(pCDH, 0, sizeof(AliRawDataHeader));
310     memset(pHLTH, 0, sizeof(AliHLTOUT::AliHLTOUTEventHeader));
311
312     if (pWriter) {
313       // copy payload
314       pWriter->Copy(&fBuffer[sizeof(AliRawDataHeader)+sizeof(AliHLTOUT::AliHLTOUTEventHeader)], 0, 0, 0, 0);
315       pHLTH->fLength=pWriter->GetTotalMemorySize();
316       // set status bit to indicate HLT payload
317       pCDH->fStatusMiniEventID|=0x1<<(AliHLTOUT::kCDHStatusFlagsOffset+AliHLTOUT::kCDHFlagsHLTPayload);
318     }
319     pHLTH->fLength+=sizeof(AliHLTOUT::AliHLTOUTEventHeader);
320     pHLTH->fEventID=eventNo;
321     // version does not really matter since we do not add decision data
322     pHLTH->fVersion=AliHLTOUT::kVersion1;
323
324     pCDH->fSize=sizeof(AliRawDataHeader)+pHLTH->fLength;
325     pCDH->fStatusMiniEventID|=0x1<<(AliHLTOUT::kCDHStatusFlagsOffset + AliHLTOUT::kCDHFlagsHLTPayload);
326     
327     pBuffer=&fBuffer[0];
328     iResult=(int)bufferSize;
329   } else {
330     pBuffer=NULL;
331     iResult=-ENOMEM;
332   }
333
334   return iResult;
335 }
336
337 int AliHLTOUTComponent::WriteDigits(int /*eventNo*/, AliRunLoader* /*runLoader*/, int /*hltddl*/, const AliHLTUInt8_t* /*pBuffer*/, unsigned int /*bufferSize*/)
338 {
339   // see header file for class documentation
340   int iResult=0;
341   return iResult;
342 }
343
344 int AliHLTOUTComponent::WriteRawFile(int eventNo, AliRunLoader* /*runLoader*/, int hltddl, const AliHLTUInt8_t* pBuffer, unsigned int bufferSize)
345 {
346   // see header file for class documentation
347   int iResult=0;
348   const char* fileName=AliDAQ::DdlFileName("HLT", hltddl);
349   assert(fileName!=NULL);
350   TString filePath;
351   filePath.Form("raw%d/", eventNo);
352   filePath+=fileName;
353   if (fileName) {
354     ios::openmode filemode=(ios::openmode)0;
355     ofstream rawfile(filePath.Data(), filemode);
356     if (rawfile.good()) {
357       if (pBuffer && bufferSize>0) {
358         rawfile.write(reinterpret_cast<const char*>(pBuffer), bufferSize);
359       } else {
360         HLTWarning("writing zero length raw data file %s");
361       }
362       HLTDebug("wrote %d byte(s) to file %s", bufferSize, filePath.Data());
363     } else {
364       HLTError("can not open file %s for writing", filePath.Data());
365       iResult=-EBADF;
366     }
367     rawfile.close();
368   }
369   return iResult;
370 }
371
372 void AliHLTOUTComponent::SetGlobalOption(unsigned int options)
373 {
374   // see header file for class documentation
375   fgOptions|=options;
376 }
377
378 void AliHLTOUTComponent::ClearGlobalOption(unsigned int options)
379 {
380   // see header file for class documentation
381   fgOptions&=~options;
382 }