3 /**************************************************************************
4 * This file is property of and copyright by the ALICE HLT Project *
5 * ALICE Experiment at CERN, All rights reserved. *
7 * Primary Authors: Matthias Richter <Matthias.Richter@ift.uib.no> *
8 * for The ALICE HLT Project. *
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 **************************************************************************/
19 /** @file AliHLTOUTComponent.cxx
20 @author Matthias Richter
22 @brief The HLTOUT data sink component similar to HLTOUT nodes
25 // see header file for class documentation
27 // refer to README to build package
29 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
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
46 /** ROOT macro for the implementation of ROOT specific class methods */
47 ClassImp(AliHLTOUTComponent)
49 AliHLTOUTComponent::AliHLTOUTComponent()
51 AliHLTOfflineDataSink(),
54 fIdFirstDDL(7680), // 0x1e<<8
58 // see header file for class documentation
60 // refer to README to build package
62 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
64 // I guess DDL definitions should never change any more
65 assert(fNofDDLs==AliDAQ::NumberOfDdls("HLT"));
66 fNofDDLs=AliDAQ::NumberOfDdls("HLT");
68 /* AliDAQ::DdlIDOffset returns wrong offset for HLT links
69 assert(fIdFirstDDL==AliDAQ::DdlIDOffset("HLT"));
70 fIdFirstDDL=AliDAQ::DdlIDOffset("HLT");
74 int AliHLTOUTComponent::fgOptions=kWriteRawFiles|kWriteDigits;
76 AliHLTOUTComponent::~AliHLTOUTComponent()
78 // see header file for class documentation
79 if (fpLibManager) delete fpLibManager;
83 const char* AliHLTOUTComponent::GetComponentID()
85 // see header file for class documentation
89 void AliHLTOUTComponent::GetInputDataTypes( vector<AliHLTComponentDataType>& list)
91 // see header file for class documentation
93 list.push_back(kAliHLTAnyDataType);
96 AliHLTComponent* AliHLTOUTComponent::Spawn()
98 // see header file for class documentation
99 return new AliHLTOUTComponent;
102 int AliHLTOUTComponent::DoInit( int argc, const char** argv )
104 // see header file for class documentation
108 for (int i=0; i<argc && iResult>=0; i++) {
110 if (argument.IsNull()) continue;
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();
120 HLTError("wrong parameter for argument %s, number expected", argument.Data());
124 HLTError("unknown argument %s", argument.Data());
130 HLTError("missing parameter for argument %s", argument.Data());
136 fpLibManager=new AliHLTHOMERLibManager;
139 for (writerNo=0; writerNo<fNofDDLs; writerNo++) {
140 AliHLTMonitoringWriter* pWriter=fpLibManager->OpenWriter();
142 HLTDebug("HOMER writer %p added", pWriter);
143 fWriters.push_back(pWriter);
145 HLTError("can nor open HOMER writer");
157 int AliHLTOUTComponent::DoDeinit()
159 // see header file for class documentation
163 AliHLTMonitoringWriterPVector::iterator element=fWriters.begin();
164 while (element!= fWriters.end()) {
166 // wanted to have a dynamic_cast<AliHLTHOMERWriter*> here, but this results into
167 // undefined symbol when loading the library
169 if (*element!=NULL) fpLibManager->DeleteWriter((AliHLTHOMERWriter*)(*element));
170 element=fWriters.erase(element);
177 int AliHLTOUTComponent::DumpEvent( const AliHLTComponentEventData& evtData,
178 const AliHLTComponentBlockData* blocks,
179 AliHLTComponentTriggerData& /*trigData*/ )
181 // see header file for class documentation
183 HLTInfo("write %d output block(s)", evtData.fBlockCnt);
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
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
211 fWriters[writerNo]->AddBlock(homerHeader, blocks[n].fPtr);
218 int AliHLTOUTComponent::FillESD(int eventNo, AliRunLoader* runLoader, AliESDEvent* /*esd*/)
220 // see header file for class documentation
222 if (fWriters.size()==0) return 0;
224 // search for the writer with the biggest data volume in order to allocate the
225 // output buffer of sufficient size
227 for (size_t i=0; i<fWriters.size(); i++) {
230 if (sorted.size()==0 || fWriters[i]->GetTotalMemorySize()<=fWriters[sorted[0]]->GetTotalMemorySize()) {
233 sorted.insert(sorted.begin(), i);
238 vector<int>::iterator ddlno=sorted.begin();
239 while (ddlno!=sorted.end()) {
240 const AliHLTUInt8_t* pBuffer=NULL;
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);
247 fWriters[*ddlno]->Clear();
253 int AliHLTOUTComponent::ShuffleWriters(AliHLTMonitoringWriterPVector &list, AliHLTUInt32_t /*size*/)
255 // see header file for class documentation
257 assert(list.size()>0);
258 if (list.size()==0) return iResult;
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())
269 if (writers.size()>0) {
271 if (writers.size()>0) {
272 // shuffle among the empty writers
275 rand.SetSeed(dt.Get()*(iResult+1));
276 i=rand.Integer(writers.size()-1);
277 assert(i>0 && i<writers.size()-1);
281 // take the writer with the least data volume
287 int AliHLTOUTComponent::FillOutputBuffer(int eventNo, AliHLTMonitoringWriter* pWriter, const AliHLTUInt8_t* &pBuffer)
289 // see header file for class documentation
291 unsigned int bufferSize=0;
293 // space for common data header
294 bufferSize+=sizeof(AliRawDataHeader);
295 assert(sizeof(AliRawDataHeader)==32);
297 // space for HLT event header
298 bufferSize+=sizeof(AliHLTOUT::AliHLTOUTEventHeader);
300 // space for payload from the writer
301 if (pWriter) bufferSize+=pWriter->GetTotalMemorySize();
303 if (bufferSize>fBuffer.size())
304 fBuffer.resize(bufferSize);
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));
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);
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;
324 pCDH->fSize=sizeof(AliRawDataHeader)+pHLTH->fLength;
325 pCDH->fStatusMiniEventID|=0x1<<(AliHLTOUT::kCDHStatusFlagsOffset + AliHLTOUT::kCDHFlagsHLTPayload);
328 iResult=(int)bufferSize;
337 int AliHLTOUTComponent::WriteDigits(int /*eventNo*/, AliRunLoader* /*runLoader*/, int /*hltddl*/, const AliHLTUInt8_t* /*pBuffer*/, unsigned int /*bufferSize*/)
339 // see header file for class documentation
344 int AliHLTOUTComponent::WriteRawFile(int eventNo, AliRunLoader* /*runLoader*/, int hltddl, const AliHLTUInt8_t* pBuffer, unsigned int bufferSize)
346 // see header file for class documentation
348 const char* fileName=AliDAQ::DdlFileName("HLT", hltddl);
349 assert(fileName!=NULL);
351 filePath.Form("raw%d/", eventNo);
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);
360 HLTWarning("writing zero length raw data file %s");
362 HLTDebug("wrote %d byte(s) to file %s", bufferSize, filePath.Data());
364 HLTError("can not open file %s for writing", filePath.Data());
372 void AliHLTOUTComponent::SetGlobalOption(unsigned int options)
374 // see header file for class documentation
378 void AliHLTOUTComponent::ClearGlobalOption(unsigned int options)
380 // see header file for class documentation