]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RAW/AliMDC.cxx
A raw-data file and Tree size limit is increased to 20GB. Wherever is needed the...
[u/mrichter/AliRoot.git] / RAW / AliMDC.cxx
1 // @(#)alimdc:$Name:  $:$Id$
2 // Author: Fons Rademakers  26/11/99
3 // Updated: Dario Favretto  15/04/2003
4
5 /**************************************************************************
6  * Copyright(c) 1998-2003, ALICE Experiment at CERN, All rights reserved. *
7  *                                                                        *
8  * Author: The ALICE Off-line Project.                                    *
9  * Contributors are mentioned in the code where appropriate.              *
10  *                                                                        *
11  * Permission to use, copy, modify and distribute this software and its   *
12  * documentation strictly for non-commercial purposes is hereby granted   *
13  * without fee, provided that the above copyright notice appears in all   *
14  * copies and that both the copyright notice and this permission notice   *
15  * appear in the supporting documentation. The authors make no claims     *
16  * about the suitability of this software for any purpose. It is          *
17  * provided "as is" without express or implied warranty.                  *
18  **************************************************************************/
19
20 /* $Id$ */
21
22 //////////////////////////////////////////////////////////////////////////
23 //                                                                      //
24 // AliMDC                                                               //
25 //                                                                      //
26 // Set of classes defining the ALICE RAW event format. The AliRawEvent  //
27 // class defines a RAW event. It consists of an AliEventHeader object   //
28 // an AliEquipmentHeader object, an AliRawData object and an array of   //
29 // sub-events, themselves also being AliRawEvents. The number of        //
30 // sub-events depends on the number of DATE LDC's.                      //
31 // The AliRawEvent objects are written to a ROOT file using different   //
32 // technologies, i.e. to local disk via AliRawDB or via rfiod using     //
33 // AliRawRFIODB or via rootd using AliRawRootdDB or to CASTOR via       //
34 // rootd using AliRawCastorDB (and for performance testing there is     //
35 // also AliRawNullDB).                                                  //
36 // The AliStats class provides statics information that is added as     //
37 // a single keyed object to each raw file.                              //
38 // The AliTagDB provides an interface to a TAG database.                //
39 // The AliMDC class is usid by the "alimdc" stand-alone program         //
40 // that reads data directly from DATE.                                  //
41 //                                                                      //
42 //////////////////////////////////////////////////////////////////////////
43
44 #include <sys/types.h>
45 #include <sys/stat.h>
46
47 #include <errno.h>
48
49 #include <TSystem.h>
50 #include <TROOT.h>
51 #include <TStopwatch.h>
52 #include <TPluginManager.h>
53
54 #include <sys/uio.h>
55 #ifdef USE_EB
56 #include "libDateEb.h"
57 #endif
58
59 #include "AliMDC.h"
60
61 #include <AliLog.h>
62 #include <AliESDEvent.h>
63
64 #include "AliRawEvent.h"
65 #include "AliRawEventHeaderBase.h"
66 #include "AliRawEquipment.h"
67 #include "AliRawEquipmentHeader.h"
68 #include "AliRawData.h"
69 #include "AliStats.h"
70 #include "AliRawDB.h"
71 #include "AliRawRFIODB.h"
72 #include "AliRawCastorDB.h"
73 #include "AliRawRootdDB.h"
74 #include "AliRawNullDB.h"
75 #include "AliTagDB.h"
76 #include "AliRawEventTag.h"
77 #include "AliFilter.h"
78
79
80
81 ClassImp(AliMDC)
82
83
84 // Filter names
85 const char* const AliMDC::fgkFilterName[kNFilters] = {"AliHoughFilter"};
86
87 //______________________________________________________________________________
88 AliMDC::AliMDC(Int_t compress, Bool_t deleteFiles, EFilterMode filterMode, 
89                Double_t maxSizeTagDB, const char* fileNameTagDB,
90                const char *guidFileFolder) :
91   fEvent(new AliRawEvent),
92   fESD(NULL),
93   fStats(NULL),
94   fRawDB(NULL),
95   fTagDB(NULL),
96   fEventTag(new AliRawEventTag),
97   fCompress(compress),
98   fDeleteFiles(deleteFiles),
99   fFilterMode(filterMode),
100   fFilters(),
101   fStop(kFALSE),
102   fIsTagDBCreated(kFALSE),
103   fMaxSizeTagDB(maxSizeTagDB),
104   fFileNameTagDB(fileNameTagDB),
105   fGuidFileFolder(guidFileFolder)
106 {
107   // Create MDC processor object.
108   // compress is the file compression mode.
109   // If deleteFiles is kTRUE the raw data files will be deleted after they
110   // were closed.
111   // If the filterMode if kFilterOff no filter algorithm will be, if it is
112   // kFilterTransparent the algorthims will be run but no events will be
113   // rejected, if it is kFilterOn the filters will be run and the event will
114   // be rejected if all filters return kFALSE.
115   // If maxSizeTagDB is greater than 0 it determines the maximal size of the
116   // tag DB and then fileNameTagDB is the directory name for the tag DB.
117   // Otherwise fileNameTagDB is the file name of the tag DB. If it is NULL
118   // no tag DB will be created.
119   // Optional 'guidFileFolder' specifies the folder in which *.guid files
120   // will be stored. In case this option is not given, the *.guid files
121   // will be written to the same folder as the raw-data files.
122
123   // Set the maximum tree size to 19GB
124   // in order to allow big raw data files
125   TTree::SetMaxTreeSize(20000000000LL);
126  
127   // This line is needed in case of a stand-alone application w/o
128   // $ROOTSYS/etc/system.rootrc file
129   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
130                                         "*",
131                                         "TStreamerInfo",
132                                         "RIO",
133                                         "TStreamerInfo()");
134
135   if (fFilterMode != kFilterOff) {
136     fESD = new AliESDEvent();
137   }
138
139   // Create the guid files folder if it does not exist
140   if (!fGuidFileFolder.IsNull()) {
141     gSystem->ResetErrno();
142     gSystem->MakeDirectory(fGuidFileFolder.Data());
143     if (gSystem->GetErrno() && gSystem->GetErrno() != EEXIST) {
144       SysError("AliMDC", "mkdir %s", fGuidFileFolder.Data());
145     }
146   }
147
148   // install SIGUSR1 handler to allow clean interrupts
149   gSystem->AddSignalHandler(new AliMDCInterruptHandler(this));
150
151   // create the high level filters
152   if (fFilterMode != kFilterOff) {
153     for (Int_t iFilter = 0; iFilter < kNFilters; iFilter++) {
154       TClass* filterClass = gROOT->GetClass(fgkFilterName[iFilter]);
155       if (!filterClass) {
156         Warning("AliMDC", "no filter class %s found", fgkFilterName[iFilter]);
157         continue;
158       }
159       AliFilter* filter = (AliFilter*) filterClass->New();
160       if (!filter) {
161         Warning("AliMDC", "creation of filter %s failed", fgkFilterName[iFilter]);
162         continue;
163       }
164       fFilters.Add(filter);
165     }
166   }
167 }
168
169 //______________________________________________________________________________
170 AliMDC::~AliMDC()
171 {
172 // destructor
173
174   fFilters.Delete();
175   if(fTagDB) delete fTagDB;
176   delete fRawDB;
177   delete fStats;
178   delete fESD;
179   delete fEvent;
180   delete fEventTag;
181 }
182  
183 //______________________________________________________________________________
184 Int_t AliMDC::Open(EWriteMode mode, const char* fileName,
185                    Double_t maxFileSize,
186                    const char* fs1, const char* fs2)
187 {
188 // open a new raw DB file
189
190   if (mode == kRFIO)
191     fRawDB = new AliRawRFIODB(fEvent, fESD, fCompress, fileName);
192   else if (mode == kROOTD)
193     fRawDB = new AliRawRootdDB(fEvent, fESD, fCompress, fileName);
194   else if (mode == kCASTOR)
195     fRawDB = new AliRawCastorDB(fEvent, fESD, fCompress, fileName);
196   else if (mode == kDEVNULL)
197     fRawDB = new AliRawNullDB(fEvent, fESD, fCompress, fileName);
198   else
199     fRawDB = new AliRawDB(fEvent, fESD, fCompress, fileName);
200   fRawDB->SetDeleteFiles(fDeleteFiles);
201
202   if (fRawDB->IsZombie()) {
203     delete fRawDB;
204     fRawDB = NULL;
205     return -1;
206   }
207
208   if (fileName == NULL) {
209     fRawDB->SetMaxSize(maxFileSize);
210     fRawDB->SetFS(fs1, fs2);
211     if (!fRawDB->Create()) {
212       delete fRawDB;
213       fRawDB = NULL;
214       return -1;
215     }
216   }
217
218   if (!fRawDB->WriteGuidFile(fGuidFileFolder)) {
219     delete fRawDB;
220     fRawDB = NULL;
221     return -2;
222   }
223
224   Info("Open", "Filling raw DB %s\n", fRawDB->GetDBName());
225
226   // Create AliStats object
227   fStats = new AliStats(fRawDB->GetDBName(), fCompress, 
228                         fFilterMode != kFilterOff);
229   return 0;
230 }
231
232 //______________________________________________________________________________
233 Int_t AliMDC::ProcessEvent(void* event, Bool_t isIovecArray)
234 {
235 // Convert the DATE event to an AliRawEvent object and store it in the raw DB,
236 // optionally also run the filter.
237 // event is either a pointer to the streamlined event 
238 // or, if isIovecArray is kTRUE, a pointer to an array of iovecs with one
239 // iovec per subevent (used by the event builder).
240 // The return value is the number of written bytes or an error code
241   const Long64_t kFileSizeErrorLevel   = 19000000000LL;
242
243   Long64_t currentFileSize = GetTotalSize();
244   if(currentFileSize > kFileSizeErrorLevel) {
245     Error("ProcessEvent", "file size (%lu) exceeds the limit "
246           , currentFileSize);
247     return kErrFileSize;
248   }
249
250   Int_t status;
251   char* data = (char*) event;
252   if (isIovecArray) data = (char*) ((iovec*) event)[0].iov_base;
253
254   // Shortcut for easy header access
255   AliRawEventHeaderBase *header = fEvent->GetHeader(data);
256
257   // Read event header
258   if ((status = header->ReadHeader(data)) != (Int_t)header->GetHeadSize()) {
259     return kErrHeader;
260   }
261
262   if (AliDebugLevel() > 2) ToAliDebug(3, header->Dump(););
263
264   // Check event type and skip "Start of Run", "End of Run",
265   // "Start of Run Files" and "End of Run Files"
266   Int_t size = header->GetEventSize() - header->GetHeadSize();
267   
268   AliDebug(1, Form("Processing %s (%d bytes)", header->GetTypeName(), size));
269
270   // Amount of data left to read for this event
271   Int_t toRead = size;
272
273   // StartOfRun, EndOfRun etc. events have no payload
274   // Nevertheless, store the event headers in the tree
275   if (toRead > 0) {
276
277     // If there is less data for this event than the next sub-event
278     // header, something is wrong. Skip to next event...
279     if (toRead < (Int_t)header->GetHeadSize()) {
280       Error("ProcessEvent", "header size (%d) exceeds number of bytes "
281             "to read (%d)", header->GetHeadSize(), toRead);
282       if (AliDebugLevel() > 0) ToAliDebug(1, header->Dump(););
283       return kErrHeaderSize;
284     }
285   
286     // Loop over all sub-events... (LDCs)
287     Int_t nsub = 1;
288     while (toRead > 0) {
289       if (isIovecArray) data = (char*) ((iovec*) event)[nsub].iov_base;
290
291       AliDebug(1, Form("reading LDC %d", nsub));
292
293       AliRawEvent *subEvent = fEvent->NextSubEvent();
294
295       // Read sub-event header
296       AliRawEventHeaderBase *subHeader = subEvent->GetHeader(data);
297       if ((status = subHeader->ReadHeader(data)) != (Int_t)subHeader->GetHeadSize()) {
298         return kErrSubHeader;
299       }
300
301       if (AliDebugLevel() > 2) ToAliDebug(3, subHeader->Dump(););
302
303       toRead -= subHeader->GetHeadSize();
304
305       Int_t rawSize = subHeader->GetEventSize() - subHeader->GetHeadSize();
306
307       // Make sure raw data less than left over bytes for current event
308       if (rawSize > toRead) {
309         Warning("ProcessEvent", "raw data size (%d) exceeds number of "
310                 "bytes to read (%d)\n", rawSize, toRead);
311         if (AliDebugLevel() > 0) ToAliDebug(1, subHeader->Dump(););
312         return kErrDataSize;
313       }
314
315       // Read Equipment Headers (in case of physics or calibration event)
316       if (header->Get("Type") == AliRawEventHeaderBase::kPhysicsEvent ||
317           header->Get("Type") == AliRawEventHeaderBase::kCalibrationEvent ||
318           header->Get("Type") == AliRawEventHeaderBase::kSystemSoftwareTriggerEvent ||
319           header->Get("Type") == AliRawEventHeaderBase::kDetectorSoftwareTriggerEvent) {
320         while (rawSize > 0) {
321           AliRawEquipment &equipment = *subEvent->NextEquipment();
322           AliRawEquipmentHeader &equipmentHeader = 
323             *equipment.GetEquipmentHeader();
324           Int_t equipHeaderSize = equipmentHeader.HeaderSize();
325           if ((status = ReadEquipmentHeader(equipmentHeader, header->DataIsSwapped(),
326                                             data)) != equipHeaderSize) {
327             return kErrEquipmentHeader;
328           }
329
330           if (AliDebugLevel() > 2) ToAliDebug(3, equipmentHeader.Dump(););
331
332           toRead  -= equipHeaderSize;
333           rawSize -= equipHeaderSize;
334
335           // Read equipment raw data
336           AliRawData &subRaw = *equipment.GetRawData();
337
338           Int_t eqSize = equipmentHeader.GetEquipmentSize() - equipHeaderSize;
339           if ((status = ReadRawData(subRaw, eqSize, data)) != eqSize) {
340             return kErrEquipment;
341           }
342           toRead  -= eqSize;
343           rawSize -= eqSize;
344
345         }
346
347       } else {  // Read only raw data but no equipment header
348         if (rawSize) {
349           AliRawEquipment &equipment = *subEvent->NextEquipment();
350           AliRawData &subRaw = *equipment.GetRawData();
351           if ((status = ReadRawData(subRaw, rawSize, data)) != rawSize) {
352             return kErrEquipment;
353           }
354           toRead  -= rawSize;
355         }
356       }
357
358       nsub++;
359     }
360   }
361
362   // High Level Event Filter
363   if (fFilterMode != kFilterOff) {
364     if (header->Get("Type") == AliRawEventHeaderBase::kPhysicsEvent ||
365         header->Get("Type") == AliRawEventHeaderBase::kCalibrationEvent ||
366         header->Get("Type") == AliRawEventHeaderBase::kSystemSoftwareTriggerEvent ||
367         header->Get("Type") == AliRawEventHeaderBase::kDetectorSoftwareTriggerEvent) {
368       Bool_t result = kFALSE;
369       for (Int_t iFilter = 0; iFilter < fFilters.GetEntriesFast(); iFilter++) {
370         AliFilter* filter = (AliFilter*) fFilters[iFilter];
371         if (!filter) continue;
372         if (filter->Filter(fEvent, fESD)) result = kTRUE;
373       }
374       if ((fFilterMode == kFilterOn) && !result) return kFilterReject;
375     }
376   }
377
378   // Set stat info for first event of this file
379   if (fRawDB->GetEvents() == 0)
380     fStats->SetFirstId(header->Get("RunNb"), header->GetP("Id")[0]);
381
382   // Store raw event in tree
383   Int_t nBytes = fRawDB->Fill();
384
385   // Fill the event tag object
386   fEventTag->SetHeader(fEvent->GetHeader());
387   fEventTag->SetGUID(fRawDB->GetDB()->GetUUID().AsString());
388   fEventTag->SetEventNumber(fRawDB->GetEvents()-1);
389
390   // Create Tag DB here only after the raw data header
391   // version was already identified
392   if (!fIsTagDBCreated) {
393     if (!fFileNameTagDB.IsNull()) {
394       if (fMaxSizeTagDB > 0) {
395         fTagDB = new AliTagDB(fEventTag, NULL);
396         fTagDB->SetMaxSize(fMaxSizeTagDB);
397         fTagDB->SetFS(fFileNameTagDB.Data());
398         if (!fTagDB->Create()) return kErrTagFile;
399       } else {
400         fTagDB = new AliTagDB(fEventTag, fFileNameTagDB.Data());
401         if (fTagDB->IsZombie()) return kErrTagFile;
402       }
403     }
404     fIsTagDBCreated = kTRUE;
405   }
406
407   // Store event tag in tree
408   if (fTagDB) fTagDB->Fill();
409
410   // Make top event object ready for next event data
411   fEvent->Reset();
412   // Clean up HLT ESD for the next event
413   if (fESD) fESD->Reset();
414
415   if(nBytes >= 0)
416     return nBytes;
417   else
418     return kErrWriting;
419 }
420
421 //______________________________________________________________________________
422 Long64_t AliMDC::GetTotalSize()
423 {
424 // return the total current raw DB file size
425
426   if (!fRawDB) return -1;
427
428   return fRawDB->GetTotalSize();
429 }
430
431 //______________________________________________________________________________
432 Int_t AliMDC::Close()
433 {
434 // close the current raw DB file
435
436   if (!fRawDB) return -1;
437
438   fRawDB->WriteStats(fStats);
439   Int_t filesize = fRawDB->Close();
440   delete fRawDB;
441   fRawDB = NULL;
442   delete fStats;
443   fStats = NULL;
444   return filesize;
445 }
446
447 //______________________________________________________________________________
448 Int_t AliMDC::Run(const char* inputFile, Bool_t loop,
449                   EWriteMode mode, Double_t maxFileSize, 
450                   const char* fs1, const char* fs2)
451 {
452   // Run the MDC processor. Read from the input stream and only return
453   // when the input gave and EOF or a fatal error occured. On success 0
454   // is returned, 1 in case of a fatality.
455   // inputFile is the name of the DATE input file; if NULL the input will
456   // be taken from the event builder.
457   // If loop is set the same input file will be reused in an infinite loop.
458   // mode specifies the type of the raw DB.
459   // maxFileSize is the maximal size of the raw DB.
460   // fs1 and fs2 are the file system locations of the raw DB.
461
462   Info("Run", "input = %s, rawdb size = %f, filter = %s, "
463        "looping = %s, compression = %d, delete files = %s",
464        inputFile ? inputFile : "event builder", maxFileSize,
465        fFilterMode == kFilterOff ? "off" : 
466        (fFilterMode == kFilterOn ? "on" : "transparent"), 
467        loop ? "yes" : "no", fCompress, fDeleteFiles ? "yes" : "no");
468
469   // Open the input file
470   Int_t fd = -1;
471   if (inputFile) {
472     if ((fd = open(inputFile, O_RDONLY)) == -1) {
473       Error("Run", "cannot open input file %s", inputFile);
474       return 1;
475     }
476   }
477
478   // Used for statistics
479   TStopwatch timer;
480   timer.Start();
481   Double_t told = 0, tnew = 0;
482   Float_t  chunkSize = maxFileSize/100, nextChunk = chunkSize;
483
484   // Create new raw DB.
485   if (fRawDB) Close();
486
487   if (Open(mode,NULL,maxFileSize,fs1,fs2) < 0) return 1;
488
489   // Process input stream
490 #ifdef USE_EB
491   Int_t eorFlag = 0;
492 #endif
493   char* event = NULL;
494   UInt_t eventSize = 0;
495   Int_t numEvents = 0;
496
497   AliRawEventHeaderBase header;
498
499   while (kTRUE) {
500
501     // If we were in looping mode stop directly after a SIGUSR1 signal
502     if (fStop) {
503       Info("Run", "Stopping loop, processed %d events", numEvents);
504       break;
505     }
506
507     if (!inputFile) {  // get data from event builder
508 #ifdef USE_EB
509       if ((eorFlag = ebEor())) break;
510       if ((event = (char*)ebGetNextEvent()) == (char*)-1) {
511         Error("Run", "error getting next event (%s)", ebGetLastError());
512         break;
513       }
514       if (event == 0) {
515         // no event, sleep for 1 second and try again
516         gSystem->Sleep(1000);
517         continue;
518       }
519 #else
520       Error("Run", "AliMDC was compiled without event builder support");
521       delete fRawDB;
522       fRawDB = NULL;
523       delete fStats;
524       fStats = NULL;
525       return 1;
526 #endif
527
528     } else {  // get data from a file
529       {
530         Int_t nrecv;
531         if ((nrecv = Read(fd, header.HeaderBaseBegin(), header.HeaderBaseSize())) !=
532             header.HeaderBaseSize()) {
533           if (nrecv == 0) {  // eof
534             if (loop) {
535               ::lseek(fd, 0, SEEK_SET);
536               continue;
537             } else {
538               break;
539             }
540           } else {
541             Error("Run", "error reading base header");
542             Close();
543             delete[] event;
544             return 1;
545           }
546         }
547       }
548       char *data = (char *)header.HeaderBaseBegin();
549       AliRawEventHeaderBase *hdr = AliRawEventHeaderBase::Create(data);
550       Int_t nrecv;
551       if ((nrecv = Read(fd, hdr->HeaderBegin(), hdr->HeaderSize())) !=
552           hdr->HeaderSize()) {
553         if (nrecv == 0) {  // eof
554           if (loop) {
555             ::lseek(fd, 0, SEEK_SET);
556             delete hdr;
557             continue;
558           } else {
559             delete hdr;
560             break;
561           }
562         } else {
563           Error("Run", "error reading header");
564           Close();
565           delete[] event;
566           delete hdr;
567           return 1;
568         }
569       }
570       if (eventSize < hdr->GetEventSize()) {
571         delete[] event;
572         eventSize = 2 * hdr->GetEventSize();
573         event = new char[eventSize];
574       }
575       memcpy(event, header.HeaderBaseBegin(), header.HeaderBaseSize());
576       memcpy(event+hdr->HeaderBaseSize(), hdr->HeaderBegin(), hdr->HeaderSize());
577       if (hdr->GetExtendedDataSize() != 0)
578         memcpy(event+hdr->HeaderBaseSize()+hdr->HeaderSize(),
579                hdr->GetExtendedData(), hdr->GetExtendedDataSize());
580       Int_t size = hdr->GetEventSize() - hdr->GetHeadSize();
581       if (Read(fd, event + hdr->GetHeadSize(), size) != size) {
582         Error("Run", "error reading data");
583         Close();
584         delete[] event;
585         delete hdr;
586         return 1;
587       }
588       delete hdr;
589     }
590
591     Int_t result = ProcessEvent(event, !inputFile);
592     if(result < -1)
593       Error("Run", "error writing data. Error code: %d",result);
594
595     if (result >= 0) {
596       numEvents++;
597       if (!(numEvents%10))
598         printf("Processed event %d (%d)\n", numEvents, fRawDB->GetEvents());
599     }
600
601     if (result > 0) {
602       // Filling time statistics
603       if (fRawDB->GetBytesWritten() > nextChunk) {
604         tnew = timer.RealTime();
605         fStats->Fill(tnew-told);
606         told = tnew;
607         timer.Continue();
608         nextChunk += chunkSize;
609       }
610
611       // Check size of raw db. If bigger than maxFileSize, close file
612       // and continue with new file.
613       if (fRawDB->GetBytesWritten() > maxFileSize) {
614
615         printf("Written raw DB at a rate of %.1f MB/s\n",
616                fRawDB->GetBytesWritten() / timer.RealTime() / 1000000.);
617
618         // Write stats object to raw db, run db, MySQL and AliEn
619         fRawDB->WriteStats(fStats);
620         delete fStats;
621         fStats = NULL;
622
623         if (!fRawDB->NextFile()) {
624           Error("Run", "error opening next raw data file");
625           Close();
626           if (inputFile) delete[] event;
627           return 1;
628         }
629
630         printf("Filling raw DB %s\n", fRawDB->GetDBName());
631         fStats = new AliStats(fRawDB->GetDBName(), fCompress, 
632                               fFilterMode != kFilterOff);
633
634         timer.Start();
635         told = 0, tnew = 0;
636         nextChunk = chunkSize;
637       }
638
639       // Check size of tag db
640       if (fTagDB && fTagDB->FileFull()) {
641         if (!fTagDB->NextFile()) {
642           delete fTagDB;
643           fTagDB = 0;
644         } else {
645           printf("Filling tag DB %s\n", fTagDB->GetDBName());
646         }
647       }
648     }
649
650     // Make top event object ready for next event data
651     //printf("Event %d has %d sub-events\n", numEvents, fEvent->GetNSubEvents());
652     //    fEvent->Reset();
653     // Clean up HLT ESD for the next event
654     //    if (fESD) fESD->Reset();
655
656     if (!inputFile) {
657 #ifdef USE_EB
658       if (!ebReleaseEvent((iovec*)event)) {
659         Error("Run", "problem releasing event (%s)", ebGetLastError());
660         break;
661       }
662 #endif
663     }
664   }
665
666   printf("Written raw DB at a rate of %.1f MB/s\n",
667          fRawDB->GetBytesWritten() / timer.RealTime() / 1000000.);
668
669   // Write stats to raw db and run db and delete stats object
670   Close();
671
672   if (!inputFile) {
673 #ifdef USE_EB
674     // Print eor flag
675     if (eorFlag) {
676       Info("Run", "event builder reported end of run (%d)", eorFlag);
677     }
678 #endif
679   } else {
680     // Close input source
681     close(fd);
682     delete [] event;
683   }
684
685   return 0;
686 }
687
688 //______________________________________________________________________________
689 Int_t AliMDC::Read(Int_t fd, void *buffer, Int_t length)
690 {
691    // Read exactly length bytes into buffer. Returns number of bytes
692    // received, returns -1 in case of error and 0 for EOF.
693
694    errno = 0;
695
696    if (fd < 0) return -1;
697
698    Int_t n, nrecv = 0;
699    char *buf = (char *)buffer;
700
701    for (n = 0; n < length; n += nrecv) {
702       if ((nrecv = read(fd, buf+n, length-n)) <= 0) {
703          if (nrecv == 0)
704             break;        // EOF
705          if (errno != EINTR)
706             SysError("Read", "read");
707          return -1;
708       }
709    }
710    return n;
711 }
712
713 //______________________________________________________________________________
714 Int_t AliMDC::ReadEquipmentHeader(AliRawEquipmentHeader &header,
715                                   Bool_t isSwapped, char*& data)
716 {
717   // Read equipment header info from DATE data stream. Returns bytes read
718   // (i.e. AliRawEquipmentHeader::HeaderSize()), -1 in case of error and
719   // 0 for EOF. If isSwapped is kTRUE the event data is byte swapped
720   // and we will swap the header to host format.
721
722   memcpy(header.HeaderBegin(), data, header.HeaderSize());
723   data += header.HeaderSize();
724
725   // Swap equipment header data if needed
726   if (isSwapped)
727     header.Swap();
728
729   if (header.GetEquipmentSize() < (UInt_t)header.HeaderSize()) {
730     Error("ReadEquipmentHeader", "invalid equipment header size");
731     // try recovery... how?
732     return -1;
733   }
734
735   return header.HeaderSize();
736 }
737
738 //______________________________________________________________________________
739 Int_t AliMDC::ReadRawData(AliRawData &raw, Int_t size, char*& data)
740 {
741   // Read raw data from DATE data stream. Returns bytes read (i.e.
742   // size), -1 in case of error and 0 for EOF.
743
744   raw.SetBuffer(data, size);
745   data += size;
746
747   return size;
748 }
749
750 //______________________________________________________________________________
751 void AliMDC::Stop()
752 {
753   // Stop the event loop
754  
755   fStop = kTRUE; 
756   if (fRawDB) fRawDB->Stop(); 
757 }
758
759