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