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