]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RAW/AliMDC.cxx
Removal of the last sprintf. Using TString instead
[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   AliDebug(1,Form("current file size is %lld bytes",currentFileSize));
245   if(currentFileSize > kFileSizeErrorLevel) {
246     Error("ProcessEvent", "file size (%lu) exceeds the limit "
247           , currentFileSize);
248     return kErrFileSize;
249   }
250
251   Int_t status;
252   char* data = (char*) event;
253   if (isIovecArray) data = (char*) ((iovec*) event)[0].iov_base;
254
255   // Shortcut for easy header access
256   AliRawEventHeaderBase *header = fEvent->GetHeader(data);
257
258   // Read event header
259   if ((status = header->ReadHeader(data)) != (Int_t)header->GetHeadSize()) {
260     return kErrHeader;
261   }
262
263   if (AliDebugLevel() > 2) ToAliDebug(3, header->Dump(););
264
265   // Check event type and skip "Start of Run", "End of Run",
266   // "Start of Run Files" and "End of Run Files"
267   Int_t size = header->GetEventSize() - header->GetHeadSize();
268   
269   AliDebug(1, Form("Processing %s (%d bytes)", header->GetTypeName(), size));
270
271   // Amount of data left to read for this event
272   Int_t toRead = size;
273
274   // StartOfRun, EndOfRun etc. events have no payload
275   // Nevertheless, store the event headers in the tree
276   if (toRead > 0) {
277
278     // If there is less data for this event than the next sub-event
279     // header, something is wrong. Skip to next event...
280     if (toRead < (Int_t)header->GetHeadSize()) {
281       Error("ProcessEvent", "header size (%d) exceeds number of bytes "
282             "to read (%d)", header->GetHeadSize(), toRead);
283       if (AliDebugLevel() > 0) ToAliDebug(1, header->Dump(););
284       return kErrHeaderSize;
285     }
286   
287     // Loop over all sub-events... (LDCs)
288     Int_t nsub = 1;
289     while (toRead > 0) {
290       if (isIovecArray) data = (char*) ((iovec*) event)[nsub].iov_base;
291
292       AliDebug(1, Form("reading LDC %d", nsub));
293
294       AliRawEvent *subEvent = fEvent->NextSubEvent();
295
296       // Read sub-event header
297       AliRawEventHeaderBase *subHeader = subEvent->GetHeader(data);
298       if ((status = subHeader->ReadHeader(data)) != (Int_t)subHeader->GetHeadSize()) {
299         return kErrSubHeader;
300       }
301
302       if (AliDebugLevel() > 2) ToAliDebug(3, subHeader->Dump(););
303
304       toRead -= subHeader->GetHeadSize();
305
306       Int_t rawSize = subHeader->GetEventSize() - subHeader->GetHeadSize();
307
308       // Make sure raw data less than left over bytes for current event
309       if (rawSize > toRead) {
310         Warning("ProcessEvent", "raw data size (%d) exceeds number of "
311                 "bytes to read (%d)\n", rawSize, toRead);
312         if (AliDebugLevel() > 0) ToAliDebug(1, subHeader->Dump(););
313         return kErrDataSize;
314       }
315
316       // Read Equipment Headers (in case of physics or calibration event)
317       if (header->Get("Type") == AliRawEventHeaderBase::kPhysicsEvent ||
318           header->Get("Type") == AliRawEventHeaderBase::kCalibrationEvent ||
319           header->Get("Type") == AliRawEventHeaderBase::kSystemSoftwareTriggerEvent ||
320           header->Get("Type") == AliRawEventHeaderBase::kDetectorSoftwareTriggerEvent) {
321         while (rawSize > 0) {
322           AliRawEquipment &equipment = *subEvent->NextEquipment();
323           AliRawEquipmentHeader &equipmentHeader = 
324             *equipment.GetEquipmentHeader();
325           Int_t equipHeaderSize = equipmentHeader.HeaderSize();
326           if ((status = ReadEquipmentHeader(equipmentHeader, header->DataIsSwapped(),
327                                             data)) != equipHeaderSize) {
328             return kErrEquipmentHeader;
329           }
330
331           if (AliDebugLevel() > 2) ToAliDebug(3, equipmentHeader.Dump(););
332
333           toRead  -= equipHeaderSize;
334           rawSize -= equipHeaderSize;
335
336           // Read equipment raw data
337           AliRawData &subRaw = *equipment.GetRawData();
338
339           Int_t eqSize = equipmentHeader.GetEquipmentSize() - equipHeaderSize;
340           if ((status = ReadRawData(subRaw, eqSize, data)) != eqSize) {
341             return kErrEquipment;
342           }
343           toRead  -= eqSize;
344           rawSize -= eqSize;
345
346         }
347
348       } else {  // Read only raw data but no equipment header
349         if (rawSize) {
350           AliRawEquipment &equipment = *subEvent->NextEquipment();
351           AliRawData &subRaw = *equipment.GetRawData();
352           if ((status = ReadRawData(subRaw, rawSize, data)) != rawSize) {
353             return kErrEquipment;
354           }
355           toRead  -= rawSize;
356         }
357       }
358
359       nsub++;
360     }
361   }
362
363   // High Level Event Filter
364   if (fFilterMode != kFilterOff) {
365     if (header->Get("Type") == AliRawEventHeaderBase::kPhysicsEvent ||
366         header->Get("Type") == AliRawEventHeaderBase::kCalibrationEvent ||
367         header->Get("Type") == AliRawEventHeaderBase::kSystemSoftwareTriggerEvent ||
368         header->Get("Type") == AliRawEventHeaderBase::kDetectorSoftwareTriggerEvent) {
369       Bool_t result = kFALSE;
370       for (Int_t iFilter = 0; iFilter < fFilters.GetEntriesFast(); iFilter++) {
371         AliFilter* filter = (AliFilter*) fFilters[iFilter];
372         if (!filter) continue;
373         if (filter->Filter(fEvent, fESD)) result = kTRUE;
374       }
375       if ((fFilterMode == kFilterOn) && !result) return kFilterReject;
376     }
377   }
378
379   // Set stat info for first event of this file
380   if (fRawDB->GetEvents() == 0)
381     fStats->SetFirstId(header->Get("RunNb"), header->GetP("Id")[0]);
382
383   // Store raw event in tree
384   Int_t nBytes = fRawDB->Fill();
385
386   // Fill the event tag object
387   fEventTag->SetHeader(fEvent->GetHeader());
388   fEventTag->SetGUID(fRawDB->GetDB()->GetUUID().AsString());
389   fEventTag->SetEventNumber(fRawDB->GetEvents()-1);
390
391   // Create Tag DB here only after the raw data header
392   // version was already identified
393   if (!fIsTagDBCreated) {
394     if (!fFileNameTagDB.IsNull()) {
395       if (fMaxSizeTagDB > 0) {
396         fTagDB = new AliTagDB(fEventTag, NULL);
397         fTagDB->SetMaxSize(fMaxSizeTagDB);
398         fTagDB->SetFS(fFileNameTagDB.Data());
399         if (!fTagDB->Create()) return kErrTagFile;
400       } else {
401         fTagDB = new AliTagDB(fEventTag, fFileNameTagDB.Data());
402         if (fTagDB->IsZombie()) return kErrTagFile;
403       }
404     }
405     fIsTagDBCreated = kTRUE;
406   }
407
408   // Store event tag in tree
409   if (fTagDB) fTagDB->Fill();
410
411   // Make top event object ready for next event data
412   fEvent->Reset();
413   // Clean up HLT ESD for the next event
414   if (fESD) fESD->Reset();
415
416   if(nBytes >= 0)
417     return nBytes;
418   else
419     return kErrWriting;
420 }
421
422 //______________________________________________________________________________
423 Long64_t AliMDC::GetTotalSize()
424 {
425 // return the total current raw DB file size
426
427   if (!fRawDB) return -1;
428
429   return fRawDB->GetTotalSize();
430 }
431
432 //______________________________________________________________________________
433 Long64_t AliMDC::Close()
434 {
435 // close the current raw DB file
436
437   if (!fRawDB) return -1;
438
439   fRawDB->WriteStats(fStats);
440   Long64_t filesize = fRawDB->Close();
441   delete fRawDB;
442   fRawDB = NULL;
443   delete fStats;
444   fStats = NULL;
445   return filesize;
446 }
447
448 //______________________________________________________________________________
449 Long64_t AliMDC::AutoSave()
450 {
451   // Auto-save the raw-data
452   // and esd (if any) trees
453
454   if (!fRawDB) return -1;
455
456   return fRawDB->AutoSave();
457 }
458
459 //______________________________________________________________________________
460 Int_t AliMDC::Run(const char* inputFile, Bool_t loop,
461                   EWriteMode mode, Double_t maxFileSize, 
462                   const char* fs1, const char* fs2)
463 {
464   // Run the MDC processor. Read from the input stream and only return
465   // when the input gave and EOF or a fatal error occured. On success 0
466   // is returned, 1 in case of a fatality.
467   // inputFile is the name of the DATE input file; if NULL the input will
468   // be taken from the event builder.
469   // If loop is set the same input file will be reused in an infinite loop.
470   // mode specifies the type of the raw DB.
471   // maxFileSize is the maximal size of the raw DB.
472   // fs1 and fs2 are the file system locations of the raw DB.
473
474   Info("Run", "input = %s, rawdb size = %f, filter = %s, "
475        "looping = %s, compression = %d, delete files = %s",
476        inputFile ? inputFile : "event builder", maxFileSize,
477        fFilterMode == kFilterOff ? "off" : 
478        (fFilterMode == kFilterOn ? "on" : "transparent"), 
479        loop ? "yes" : "no", fCompress, fDeleteFiles ? "yes" : "no");
480
481   // Open the input file
482   Int_t fd = -1;
483   if (inputFile) {
484     if ((fd = open(inputFile, O_RDONLY)) == -1) {
485       Error("Run", "cannot open input file %s", inputFile);
486       return 1;
487     }
488   }
489
490   // Used for statistics
491   TStopwatch timer;
492   timer.Start();
493   Double_t told = 0, tnew = 0;
494   Float_t  chunkSize = maxFileSize/100, nextChunk = chunkSize;
495
496   // Create new raw DB.
497   if (fRawDB) Close();
498
499   if (Open(mode,NULL,maxFileSize,fs1,fs2) < 0) return 1;
500
501   // Process input stream
502 #ifdef USE_EB
503   Int_t eorFlag = 0;
504 #endif
505   char* event = NULL;
506   UInt_t eventSize = 0;
507   Int_t numEvents = 0;
508
509   AliRawEventHeaderBase header;
510
511   while (kTRUE) {
512
513     // If we were in looping mode stop directly after a SIGUSR1 signal
514     if (fStop) {
515       Info("Run", "Stopping loop, processed %d events", numEvents);
516       break;
517     }
518
519     if (!inputFile) {  // get data from event builder
520 #ifdef USE_EB
521       if ((eorFlag = ebEor())) break;
522       if ((event = (char*)ebGetNextEvent()) == (char*)-1) {
523         Error("Run", "error getting next event (%s)", ebGetLastError());
524         break;
525       }
526       if (event == 0) {
527         // no event, sleep for 1 second and try again
528         gSystem->Sleep(1000);
529         continue;
530       }
531 #else
532       Error("Run", "AliMDC was compiled without event builder support");
533       delete fRawDB;
534       fRawDB = NULL;
535       delete fStats;
536       fStats = NULL;
537       return 1;
538 #endif
539
540     } else {  // get data from a file
541       {
542         Int_t nrecv;
543         if ((nrecv = Read(fd, header.HeaderBaseBegin(), header.HeaderBaseSize())) !=
544             header.HeaderBaseSize()) {
545           if (nrecv == 0) {  // eof
546             if (loop) {
547               ::lseek(fd, 0, SEEK_SET);
548               continue;
549             } else {
550               break;
551             }
552           } else {
553             Error("Run", "error reading base header");
554             Close();
555             delete[] event;
556             return 1;
557           }
558         }
559       }
560       char *data = (char *)header.HeaderBaseBegin();
561       AliRawEventHeaderBase *hdr = AliRawEventHeaderBase::Create(data);
562       Int_t nrecv;
563       if ((nrecv = Read(fd, hdr->HeaderBegin(), hdr->HeaderSize())) !=
564           hdr->HeaderSize()) {
565         if (nrecv == 0) {  // eof
566           if (loop) {
567             ::lseek(fd, 0, SEEK_SET);
568             delete hdr;
569             continue;
570           } else {
571             delete hdr;
572             break;
573           }
574         } else {
575           Error("Run", "error reading header");
576           Close();
577           delete[] event;
578           delete hdr;
579           return 1;
580         }
581       }
582       if (eventSize < hdr->GetEventSize()) {
583         delete[] event;
584         eventSize = 2 * hdr->GetEventSize();
585         event = new char[eventSize];
586       }
587       memcpy(event, header.HeaderBaseBegin(), header.HeaderBaseSize());
588       memcpy(event+hdr->HeaderBaseSize(), hdr->HeaderBegin(), hdr->HeaderSize());
589       if (hdr->GetExtendedDataSize() != 0)
590         memcpy(event+hdr->HeaderBaseSize()+hdr->HeaderSize(),
591                hdr->GetExtendedData(), hdr->GetExtendedDataSize());
592       Int_t size = hdr->GetEventSize() - hdr->GetHeadSize();
593       if (Read(fd, event + hdr->GetHeadSize(), size) != size) {
594         Error("Run", "error reading data");
595         Close();
596         delete[] event;
597         delete hdr;
598         return 1;
599       }
600       delete hdr;
601     }
602
603     Int_t result = ProcessEvent(event, !inputFile);
604     if(result < -1)
605       Error("Run", "error writing data. Error code: %d",result);
606
607     if (result >= 0) {
608       numEvents++;
609       if (!(numEvents%10))
610         printf("Processed event %d (%d)\n", numEvents, fRawDB->GetEvents());
611     }
612
613     if (result > 0) {
614       // Filling time statistics
615       if (fRawDB->GetBytesWritten() > nextChunk) {
616         tnew = timer.RealTime();
617         fStats->Fill(tnew-told);
618         told = tnew;
619         timer.Continue();
620         nextChunk += chunkSize;
621       }
622
623       // Check size of raw db. If bigger than maxFileSize, close file
624       // and continue with new file.
625       if (fRawDB->GetBytesWritten() > maxFileSize) {
626
627         printf("Written raw DB at a rate of %.1f MB/s\n",
628                fRawDB->GetBytesWritten() / timer.RealTime() / 1000000.);
629
630         // Write stats object to raw db, run db, MySQL and AliEn
631         fRawDB->WriteStats(fStats);
632         delete fStats;
633         fStats = NULL;
634
635         if (!fRawDB->NextFile()) {
636           Error("Run", "error opening next raw data file");
637           Close();
638           if (inputFile) delete[] event;
639           return 1;
640         }
641
642         printf("Filling raw DB %s\n", fRawDB->GetDBName());
643         fStats = new AliStats(fRawDB->GetDBName(), fCompress, 
644                               fFilterMode != kFilterOff);
645
646         timer.Start();
647         told = 0, tnew = 0;
648         nextChunk = chunkSize;
649       }
650
651       // Check size of tag db
652       if (fTagDB && fTagDB->FileFull()) {
653         if (!fTagDB->NextFile()) {
654           delete fTagDB;
655           fTagDB = 0;
656         } else {
657           printf("Filling tag DB %s\n", fTagDB->GetDBName());
658         }
659       }
660     }
661
662     // Make top event object ready for next event data
663     //printf("Event %d has %d sub-events\n", numEvents, fEvent->GetNSubEvents());
664     //    fEvent->Reset();
665     // Clean up HLT ESD for the next event
666     //    if (fESD) fESD->Reset();
667
668     if (!inputFile) {
669 #ifdef USE_EB
670       if (!ebReleaseEvent((iovec*)event)) {
671         Error("Run", "problem releasing event (%s)", ebGetLastError());
672         break;
673       }
674 #endif
675     }
676   }
677
678   printf("Written raw DB at a rate of %.1f MB/s\n",
679          fRawDB->GetBytesWritten() / timer.RealTime() / 1000000.);
680
681   // Write stats to raw db and run db and delete stats object
682   Close();
683
684   if (!inputFile) {
685 #ifdef USE_EB
686     // Print eor flag
687     if (eorFlag) {
688       Info("Run", "event builder reported end of run (%d)", eorFlag);
689     }
690 #endif
691   } else {
692     // Close input source
693     close(fd);
694     delete [] event;
695   }
696
697   return 0;
698 }
699
700 //______________________________________________________________________________
701 Int_t AliMDC::Read(Int_t fd, void *buffer, Int_t length)
702 {
703    // Read exactly length bytes into buffer. Returns number of bytes
704    // received, returns -1 in case of error and 0 for EOF.
705
706    errno = 0;
707
708    if (fd < 0) return -1;
709
710    Int_t n, nrecv = 0;
711    char *buf = (char *)buffer;
712
713    for (n = 0; n < length; n += nrecv) {
714       if ((nrecv = read(fd, buf+n, length-n)) <= 0) {
715          if (nrecv == 0)
716             break;        // EOF
717          if (errno != EINTR)
718             SysError("Read", "read");
719          return -1;
720       }
721    }
722    return n;
723 }
724
725 //______________________________________________________________________________
726 Int_t AliMDC::ReadEquipmentHeader(AliRawEquipmentHeader &header,
727                                   Bool_t isSwapped, char*& data)
728 {
729   // Read equipment header info from DATE data stream. Returns bytes read
730   // (i.e. AliRawEquipmentHeader::HeaderSize()), -1 in case of error and
731   // 0 for EOF. If isSwapped is kTRUE the event data is byte swapped
732   // and we will swap the header to host format.
733
734   memcpy(header.HeaderBegin(), data, header.HeaderSize());
735   data += header.HeaderSize();
736
737   // Swap equipment header data if needed
738   if (isSwapped)
739     header.Swap();
740
741   if (header.GetEquipmentSize() < (UInt_t)header.HeaderSize()) {
742     Error("ReadEquipmentHeader", "invalid equipment header size");
743     // try recovery... how?
744     return -1;
745   }
746
747   return header.HeaderSize();
748 }
749
750 //______________________________________________________________________________
751 Int_t AliMDC::ReadRawData(AliRawData &raw, Int_t size, char*& data)
752 {
753   // Read raw data from DATE data stream. Returns bytes read (i.e.
754   // size), -1 in case of error and 0 for EOF.
755
756   raw.SetBuffer(data, size);
757   data += size;
758
759   return size;
760 }
761
762 //______________________________________________________________________________
763 void AliMDC::Stop()
764 {
765   // Stop the event loop
766  
767   fStop = kTRUE; 
768   if (fRawDB) fRawDB->Stop(); 
769 }
770
771