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