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