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