]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RAW/AliMDC.cxx
New class AliESDEvent, backward compatibility with the old AliESD (Christian)
[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   // This line is needed in case of a stand-alone application w/o
124   // $ROOTSYS/etc/system.rootrc file
125   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
126                                         "*",
127                                         "TStreamerInfo",
128                                         "RIO",
129                                         "TStreamerInfo()");
130
131   if (fFilterMode != kFilterOff) {
132     fESD = new AliESDEvent();
133   }
134
135   // Create the guid files folder if it does not exist
136   if (fGuidFileFolder) {
137     gSystem->ResetErrno();
138     gSystem->MakeDirectory(fGuidFileFolder);
139     if (gSystem->GetErrno() && gSystem->GetErrno() != EEXIST) {
140       SysError("AliMDC", "mkdir %s", fGuidFileFolder);
141     }
142   }
143
144   // install SIGUSR1 handler to allow clean interrupts
145   gSystem->AddSignalHandler(new AliMDCInterruptHandler(this));
146
147   // create the high level filters
148   if (fFilterMode != kFilterOff) {
149     for (Int_t iFilter = 0; iFilter < kNFilters; iFilter++) {
150       TClass* filterClass = gROOT->GetClass(fgkFilterName[iFilter]);
151       if (!filterClass) {
152         Warning("AliMDC", "no filter class %s found", fgkFilterName[iFilter]);
153         continue;
154       }
155       AliFilter* filter = (AliFilter*) filterClass->New();
156       if (!filter) {
157         Warning("AliMDC", "creation of filter %s failed", fgkFilterName[iFilter]);
158         continue;
159       }
160       fFilters.Add(filter);
161     }
162   }
163 }
164
165 //______________________________________________________________________________
166 AliMDC::~AliMDC()
167 {
168 // destructor
169
170   fFilters.Delete();
171   if(fTagDB) delete fTagDB;
172   delete fRawDB;
173   delete fStats;
174   delete fESD;
175   delete fEvent;
176   delete fEventTag;
177 }
178  
179 //______________________________________________________________________________
180 Int_t AliMDC::Open(EWriteMode mode, const char* fileName,
181                    Double_t maxFileSize,
182                    const char* fs1, const char* fs2)
183 {
184 // open a new raw DB file
185
186   if (mode == kRFIO)
187     fRawDB = new AliRawRFIODB(fEvent, fESD, fCompress, fileName);
188   else if (mode == kROOTD)
189     fRawDB = new AliRawRootdDB(fEvent, fESD, fCompress, fileName);
190   else if (mode == kCASTOR)
191     fRawDB = new AliRawCastorDB(fEvent, fESD, fCompress, fileName);
192   else if (mode == kDEVNULL)
193     fRawDB = new AliRawNullDB(fEvent, fESD, fCompress, fileName);
194   else
195     fRawDB = new AliRawDB(fEvent, fESD, fCompress, fileName);
196   fRawDB->SetDeleteFiles(fDeleteFiles);
197
198   if (fRawDB->IsZombie()) {
199     delete fRawDB;
200     fRawDB = NULL;
201     return -1;
202   }
203
204   if (fileName == NULL) {
205     fRawDB->SetMaxSize(maxFileSize);
206     fRawDB->SetFS(fs1, fs2);
207     if (!fRawDB->Create()) {
208       delete fRawDB;
209       fRawDB = NULL;
210       return -1;
211     }
212   }
213
214   if (!fRawDB->WriteGuidFile(fGuidFileFolder)) {
215     delete fRawDB;
216     fRawDB = NULL;
217     return -2;
218   }
219
220   Info("Open", "Filling raw DB %s\n", fRawDB->GetDBName());
221
222   // Create AliStats object
223   fStats = new AliStats(fRawDB->GetDBName(), fCompress, 
224                         fFilterMode != kFilterOff);
225   return 0;
226 }
227
228 //______________________________________________________________________________
229 Int_t AliMDC::ProcessEvent(void* event, Bool_t isIovecArray)
230 {
231 // Convert the DATE event to an AliRawEvent object and store it in the raw DB,
232 // optionally also run the filter.
233 // event is either a pointer to the streamlined event 
234 // or, if isIovecArray is kTRUE, a pointer to an array of iovecs with one
235 // iovec per subevent (used by the event builder).
236 // The return value is the number of written bytes or an error code
237   const UInt_t kFileSizeErrorLevel   = 1900000000;
238
239   UInt_t currentFileSize = GetTotalSize();
240   if(currentFileSize > kFileSizeErrorLevel) {
241     Error("ProcessEvent", "file size (%u) exceeds the limit "
242           , currentFileSize);
243     return kErrFileSize;
244   }
245
246   Int_t status;
247   char* data = (char*) event;
248   if (isIovecArray) data = (char*) ((iovec*) event)[0].iov_base;
249
250   // Shortcut for easy header access
251   AliRawEventHeaderBase *header = fEvent->GetHeader(data);
252
253   // Read event header
254   if ((status = header->ReadHeader(data)) != (Int_t)header->GetHeadSize()) {
255     return kErrHeader;
256   }
257
258   if (AliDebugLevel() > 2) ToAliDebug(3, header->Dump(););
259
260   // Check event type and skip "Start of Run", "End of Run",
261   // "Start of Run Files" and "End of Run Files"
262   Int_t size = header->GetEventSize() - header->GetHeadSize();
263   
264   AliDebug(1, Form("Processing %s (%d bytes)", header->GetTypeName(), size));
265
266   // Amount of data left to read for this event
267   Int_t toRead = size;
268
269   // StartOfRun, EndOfRun etc. events have no payload
270   // Nevertheless, store the event headers in the tree
271   if (toRead > 0) {
272
273     // If there is less data for this event than the next sub-event
274     // header, something is wrong. Skip to next event...
275     if (toRead < (Int_t)header->GetHeadSize()) {
276       Error("ProcessEvent", "header size (%d) exceeds number of bytes "
277             "to read (%d)", header->GetHeadSize(), toRead);
278       if (AliDebugLevel() > 0) ToAliDebug(1, header->Dump(););
279       return kErrHeaderSize;
280     }
281   
282     // Loop over all sub-events... (LDCs)
283     Int_t nsub = 1;
284     while (toRead > 0) {
285       if (isIovecArray) data = (char*) ((iovec*) event)[nsub].iov_base;
286
287       AliDebug(1, Form("reading LDC %d", nsub));
288
289       AliRawEvent *subEvent = fEvent->NextSubEvent();
290
291       // Read sub-event header
292       AliRawEventHeaderBase *subHeader = subEvent->GetHeader(data);
293       if ((status = subHeader->ReadHeader(data)) != (Int_t)subHeader->GetHeadSize()) {
294         return kErrSubHeader;
295       }
296
297       if (AliDebugLevel() > 2) ToAliDebug(3, subHeader->Dump(););
298
299       toRead -= subHeader->GetHeadSize();
300
301       Int_t rawSize = subHeader->GetEventSize() - subHeader->GetHeadSize();
302
303       // Make sure raw data less than left over bytes for current event
304       if (rawSize > toRead) {
305         Warning("ProcessEvent", "raw data size (%d) exceeds number of "
306                 "bytes to read (%d)\n", rawSize, toRead);
307         if (AliDebugLevel() > 0) ToAliDebug(1, subHeader->Dump(););
308         return kErrDataSize;
309       }
310
311       // Read Equipment Headers (in case of physics or calibration event)
312       if (header->Get("Type") == AliRawEventHeaderBase::kPhysicsEvent ||
313           header->Get("Type") == AliRawEventHeaderBase::kCalibrationEvent ||
314           header->Get("Type") == AliRawEventHeaderBase::kSystemSoftwareTriggerEvent ||
315           header->Get("Type") == AliRawEventHeaderBase::kDetectorSoftwareTriggerEvent) {
316         while (rawSize > 0) {
317           AliRawEquipment &equipment = *subEvent->NextEquipment();
318           AliRawEquipmentHeader &equipmentHeader = 
319             *equipment.GetEquipmentHeader();
320           Int_t equipHeaderSize = equipmentHeader.HeaderSize();
321           if ((status = ReadEquipmentHeader(equipmentHeader, header->DataIsSwapped(),
322                                             data)) != equipHeaderSize) {
323             return kErrEquipmentHeader;
324           }
325
326           if (AliDebugLevel() > 2) ToAliDebug(3, equipmentHeader.Dump(););
327
328           toRead  -= equipHeaderSize;
329           rawSize -= equipHeaderSize;
330
331           // Read equipment raw data
332           AliRawData &subRaw = *equipment.GetRawData();
333
334           Int_t eqSize = equipmentHeader.GetEquipmentSize() - equipHeaderSize;
335           if ((status = ReadRawData(subRaw, eqSize, data)) != eqSize) {
336             return kErrEquipment;
337           }
338           toRead  -= eqSize;
339           rawSize -= eqSize;
340
341         }
342
343       } else {  // Read only raw data but no equipment header
344         if (rawSize) {
345           AliRawEquipment &equipment = *subEvent->NextEquipment();
346           AliRawData &subRaw = *equipment.GetRawData();
347           if ((status = ReadRawData(subRaw, rawSize, data)) != rawSize) {
348             return kErrEquipment;
349           }
350           toRead  -= rawSize;
351         }
352       }
353
354       nsub++;
355     }
356   }
357
358   // High Level Event Filter
359   if (fFilterMode != kFilterOff) {
360     if (header->Get("Type") == AliRawEventHeaderBase::kPhysicsEvent ||
361         header->Get("Type") == AliRawEventHeaderBase::kCalibrationEvent ||
362         header->Get("Type") == AliRawEventHeaderBase::kSystemSoftwareTriggerEvent ||
363         header->Get("Type") == AliRawEventHeaderBase::kDetectorSoftwareTriggerEvent) {
364       Bool_t result = kFALSE;
365       for (Int_t iFilter = 0; iFilter < fFilters.GetEntriesFast(); iFilter++) {
366         AliFilter* filter = (AliFilter*) fFilters[iFilter];
367         if (!filter) continue;
368         if (filter->Filter(fEvent, fESD)) result = kTRUE;
369       }
370       if ((fFilterMode == kFilterOn) && !result) return kFilterReject;
371     }
372   }
373
374   // Set stat info for first event of this file
375   if (fRawDB->GetEvents() == 0)
376     fStats->SetFirstId(header->Get("RunNb"), header->GetP("Id")[0]);
377
378   // Store raw event in tree
379   Int_t nBytes = fRawDB->Fill();
380
381   // Fill the event tag object
382   fEventTag->SetHeader(fEvent->GetHeader());
383   fEventTag->SetGUID(fRawDB->GetDB()->GetUUID().AsString());
384   fEventTag->SetEventNumber(fRawDB->GetEvents()-1);
385
386   // Create Tag DB here only after the raw data header
387   // version was already identified
388   if (!fIsTagDBCreated) {
389     if (fFileNameTagDB) {
390       if (fMaxSizeTagDB > 0) {
391         fTagDB = new AliTagDB(fEventTag, NULL);
392         fTagDB->SetMaxSize(fMaxSizeTagDB);
393         fTagDB->SetFS(fFileNameTagDB);
394         if (!fTagDB->Create()) return kErrTagFile;
395       } else {
396         fTagDB = new AliTagDB(fEventTag, fFileNameTagDB);
397         if (fTagDB->IsZombie()) return kErrTagFile;
398       }
399     }
400     fIsTagDBCreated = kTRUE;
401   }
402
403   // Store event tag in tree
404   if (fTagDB) fTagDB->Fill();
405
406   // Make top event object ready for next event data
407   fEvent->Reset();
408   // Clean up HLT ESD for the next event
409   if (fESD) fESD->Reset();
410
411   if(nBytes >= 0)
412     return nBytes;
413   else
414     return kErrWriting;
415 }
416
417 //______________________________________________________________________________
418 Int_t AliMDC::GetTotalSize()
419 {
420 // return the total current raw DB file size
421
422   if (!fRawDB) return -1;
423
424   return fRawDB->GetTotalSize();
425 }
426
427 //______________________________________________________________________________
428 Int_t AliMDC::Close()
429 {
430 // close the current raw DB file
431
432   if (!fRawDB) return -1;
433
434   fRawDB->WriteStats(fStats);
435   Int_t filesize = fRawDB->Close();
436   delete fRawDB;
437   fRawDB = NULL;
438   delete fStats;
439   fStats = NULL;
440   return filesize;
441 }
442
443 //______________________________________________________________________________
444 Int_t AliMDC::Run(const char* inputFile, Bool_t loop,
445                   EWriteMode mode, Double_t maxFileSize, 
446                   const char* fs1, const char* fs2)
447 {
448   // Run the MDC processor. Read from the input stream and only return
449   // when the input gave and EOF or a fatal error occured. On success 0
450   // is returned, 1 in case of a fatality.
451   // inputFile is the name of the DATE input file; if NULL the input will
452   // be taken from the event builder.
453   // If loop is set the same input file will be reused in an infinite loop.
454   // mode specifies the type of the raw DB.
455   // maxFileSize is the maximal size of the raw DB.
456   // fs1 and fs2 are the file system locations of the raw DB.
457
458   Info("Run", "input = %s, rawdb size = %f, filter = %s, "
459        "looping = %s, compression = %d, delete files = %s",
460        inputFile ? inputFile : "event builder", maxFileSize,
461        fFilterMode == kFilterOff ? "off" : 
462        (fFilterMode == kFilterOn ? "on" : "transparent"), 
463        loop ? "yes" : "no", fCompress, fDeleteFiles ? "yes" : "no");
464
465   // Open the input file
466   Int_t fd = -1;
467   if (inputFile) {
468     if ((fd = open(inputFile, O_RDONLY)) == -1) {
469       Error("Run", "cannot open input file %s", inputFile);
470       return 1;
471     }
472   }
473
474   // Used for statistics
475   TStopwatch timer;
476   timer.Start();
477   Double_t told = 0, tnew = 0;
478   Float_t  chunkSize = maxFileSize/100, nextChunk = chunkSize;
479
480   // Create new raw DB.
481   if (fRawDB) Close();
482
483   if (Open(mode,NULL,maxFileSize,fs1,fs2) < 0) return 1;
484
485   // Process input stream
486 #ifdef USE_EB
487   Int_t eorFlag = 0;
488 #endif
489   char* event = NULL;
490   UInt_t eventSize = 0;
491   Int_t numEvents = 0;
492
493   AliRawEventHeaderBase header;
494
495   while (kTRUE) {
496
497     // If we were in looping mode stop directly after a SIGUSR1 signal
498     if (fStop) {
499       Info("Run", "Stopping loop, processed %d events", numEvents);
500       break;
501     }
502
503     if (!inputFile) {  // get data from event builder
504 #ifdef USE_EB
505       if ((eorFlag = ebEor())) break;
506       if ((event = (char*)ebGetNextEvent()) == (char*)-1) {
507         Error("Run", "error getting next event (%s)", ebGetLastError());
508         break;
509       }
510       if (event == 0) {
511         // no event, sleep for 1 second and try again
512         gSystem->Sleep(1000);
513         continue;
514       }
515 #else
516       Error("Run", "AliMDC was compiled without event builder support");
517       delete fRawDB;
518       fRawDB = NULL;
519       delete fStats;
520       fStats = NULL;
521       return 1;
522 #endif
523
524     } else {  // get data from a file
525       {
526         Int_t nrecv;
527         if ((nrecv = Read(fd, header.HeaderBaseBegin(), header.HeaderBaseSize())) !=
528             header.HeaderBaseSize()) {
529           if (nrecv == 0) {  // eof
530             if (loop) {
531               ::lseek(fd, 0, SEEK_SET);
532               continue;
533             } else {
534               break;
535             }
536           } else {
537             Error("Run", "error reading base header");
538             Close();
539             delete[] event;
540             return 1;
541           }
542         }
543       }
544       char *data = (char *)header.HeaderBaseBegin();
545       AliRawEventHeaderBase *hdr = AliRawEventHeaderBase::Create(data);
546       Int_t nrecv;
547       if ((nrecv = Read(fd, hdr->HeaderBegin(), hdr->HeaderSize())) !=
548           hdr->HeaderSize()) {
549         if (nrecv == 0) {  // eof
550           if (loop) {
551             ::lseek(fd, 0, SEEK_SET);
552             delete hdr;
553             continue;
554           } else {
555             delete hdr;
556             break;
557           }
558         } else {
559           Error("Run", "error reading header");
560           Close();
561           delete[] event;
562           delete hdr;
563           return 1;
564         }
565       }
566       if (eventSize < hdr->GetEventSize()) {
567         delete[] event;
568         eventSize = 2 * hdr->GetEventSize();
569         event = new char[eventSize];
570       }
571       memcpy(event, hdr->HeaderBaseBegin(), hdr->HeaderBaseSize());
572       memcpy(event+hdr->HeaderBaseSize(), hdr->HeaderBegin(), hdr->HeaderSize());
573       if (hdr->GetExtendedDataSize() != 0)
574         memcpy(event+hdr->HeaderBaseSize()+hdr->HeaderSize(),
575                hdr->GetExtendedData(), hdr->GetExtendedDataSize());
576       Int_t size = hdr->GetEventSize() - hdr->GetHeadSize();
577       if (Read(fd, event + hdr->GetHeadSize(), size) != size) {
578         Error("Run", "error reading data");
579         Close();
580         delete[] event;
581         delete hdr;
582         return 1;
583       }
584       delete hdr;
585     }
586
587     Int_t result = ProcessEvent(event, !inputFile);
588     if(result < -1)
589       Error("Run", "error writing data. Error code: %d",result);
590
591     if (result >= 0) {
592       numEvents++;
593       if (!(numEvents%10))
594         printf("Processed event %d (%d)\n", numEvents, fRawDB->GetEvents());
595     }
596
597     if (result > 0) {
598       // Filling time statistics
599       if (fRawDB->GetBytesWritten() > nextChunk) {
600         tnew = timer.RealTime();
601         fStats->Fill(tnew-told);
602         told = tnew;
603         timer.Continue();
604         nextChunk += chunkSize;
605       }
606
607       // Check size of raw db. If bigger than maxFileSize, close file
608       // and continue with new file.
609       if (fRawDB->GetBytesWritten() > maxFileSize) {
610
611         printf("Written raw DB at a rate of %.1f MB/s\n",
612                fRawDB->GetBytesWritten() / timer.RealTime() / 1000000.);
613
614         // Write stats object to raw db, run db, MySQL and AliEn
615         fRawDB->WriteStats(fStats);
616         delete fStats;
617         fStats = NULL;
618
619         if (!fRawDB->NextFile()) {
620           Error("Run", "error opening next raw data file");
621           Close();
622           if (inputFile) delete[] event;
623           return 1;
624         }
625
626         printf("Filling raw DB %s\n", fRawDB->GetDBName());
627         fStats = new AliStats(fRawDB->GetDBName(), fCompress, 
628                               fFilterMode != kFilterOff);
629
630         timer.Start();
631         told = 0, tnew = 0;
632         nextChunk = chunkSize;
633       }
634
635       // Check size of tag db
636       if (fTagDB && fTagDB->FileFull()) {
637         if (!fTagDB->NextFile()) {
638           delete fTagDB;
639           fTagDB = 0;
640         } else {
641           printf("Filling tag DB %s\n", fTagDB->GetDBName());
642         }
643       }
644     }
645
646     // Make top event object ready for next event data
647     //printf("Event %d has %d sub-events\n", numEvents, fEvent->GetNSubEvents());
648     //    fEvent->Reset();
649     // Clean up HLT ESD for the next event
650     //    if (fESD) fESD->Reset();
651
652     if (!inputFile) {
653 #ifdef USE_EB
654       if (!ebReleaseEvent((iovec*)event)) {
655         Error("Run", "problem releasing event (%s)", ebGetLastError());
656         break;
657       }
658 #endif
659     }
660   }
661
662   printf("Written raw DB at a rate of %.1f MB/s\n",
663          fRawDB->GetBytesWritten() / timer.RealTime() / 1000000.);
664
665   // Write stats to raw db and run db and delete stats object
666   Close();
667
668   if (!inputFile) {
669 #ifdef USE_EB
670     // Print eor flag
671     if (eorFlag) {
672       Info("Run", "event builder reported end of run (%d)", eorFlag);
673     }
674 #endif
675   } else {
676     // Close input source
677     close(fd);
678     delete [] event;
679   }
680
681   return 0;
682 }
683
684 //______________________________________________________________________________
685 Int_t AliMDC::Read(Int_t fd, void *buffer, Int_t length)
686 {
687    // Read exactly length bytes into buffer. Returns number of bytes
688    // received, returns -1 in case of error and 0 for EOF.
689
690    errno = 0;
691
692    if (fd < 0) return -1;
693
694    Int_t n, nrecv = 0;
695    char *buf = (char *)buffer;
696
697    for (n = 0; n < length; n += nrecv) {
698       if ((nrecv = read(fd, buf+n, length-n)) <= 0) {
699          if (nrecv == 0)
700             break;        // EOF
701          if (errno != EINTR)
702             SysError("Read", "read");
703          return -1;
704       }
705    }
706    return n;
707 }
708
709 //______________________________________________________________________________
710 Int_t AliMDC::ReadEquipmentHeader(AliRawEquipmentHeader &header,
711                                   Bool_t isSwapped, char*& data)
712 {
713   // Read equipment header info from DATE data stream. Returns bytes read
714   // (i.e. AliRawEquipmentHeader::HeaderSize()), -1 in case of error and
715   // 0 for EOF. If isSwapped is kTRUE the event data is byte swapped
716   // and we will swap the header to host format.
717
718   memcpy(header.HeaderBegin(), data, header.HeaderSize());
719   data += header.HeaderSize();
720
721   // Swap equipment header data if needed
722   if (isSwapped)
723     header.Swap();
724
725   if (header.GetEquipmentSize() < (UInt_t)header.HeaderSize()) {
726     Error("ReadEquipmentHeader", "invalid equipment header size");
727     // try recovery... how?
728     return -1;
729   }
730
731   return header.HeaderSize();
732 }
733
734 //______________________________________________________________________________
735 Int_t AliMDC::ReadRawData(AliRawData &raw, Int_t size, char*& data)
736 {
737   // Read raw data from DATE data stream. Returns bytes read (i.e.
738   // size), -1 in case of error and 0 for EOF.
739
740   raw.SetBuffer(data, size);
741   data += size;
742
743   return size;
744 }
745
746 //______________________________________________________________________________
747 void AliMDC::Stop()
748 {
749   // Stop the event loop
750  
751   fStop = kTRUE; 
752   if (fRawDB) fRawDB->Stop(); 
753 }
754
755