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