]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RAW/AliRawReaderRoot.cxx
Coverity 22681
[u/mrichter/AliRoot.git] / RAW / AliRawReaderRoot.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 ///////////////////////////////////////////////////////////////////////////////
19 ///
20 /// This is a class for reading raw data from a root file.
21 ///
22 /// The root file is expected to contain a tree of name "RAW" with
23 /// a branch of name "rawevent" which contains objects of type
24 /// AliRawVEvent.
25 /// 
26 /// The file name and the event number are arguments of the constructor
27 /// of AliRawReaderRoot.
28 ///
29 ///////////////////////////////////////////////////////////////////////////////
30
31 #include <TFile.h>
32 #include <TTree.h>
33 #include <TTreeIndex.h>
34 #include "AliRawReaderRoot.h"
35 #include "AliRawVEvent.h"
36 #include "AliRawEventHeaderBase.h"
37 #include "AliRawVEquipment.h"
38 #include "AliRawEquipmentHeader.h"
39 #include "AliRawData.h"
40
41
42 ClassImp(AliRawReaderRoot)
43 Bool_t AliRawReaderRoot::fgUseOrder = kFALSE;
44
45
46 AliRawReaderRoot::AliRawReaderRoot() :
47   fFile(NULL),
48   fBranch(NULL),
49   fEventIndex(-1),
50   fEvent(NULL),
51   fEventHeader(NULL),
52   fSubEventIndex(0),
53   fSubEvent(NULL),
54   fEquipmentIndex(0),
55   fEquipment(NULL),
56   fRawData(NULL),
57   fPosition(NULL),
58   fEnd(NULL),
59   fIndex(0x0)
60 {
61 // default constructor
62
63 }
64
65 AliRawReaderRoot::AliRawReaderRoot(const char* fileName, Int_t eventNumber) :
66   fFile(NULL),
67   fBranch(NULL),
68   fEventIndex(eventNumber),
69   fEvent(NULL),
70   fEventHeader(NULL),
71   fSubEventIndex(0),
72   fSubEvent(NULL),
73   fEquipmentIndex(0),
74   fEquipment(NULL),
75   fRawData(NULL),
76   fPosition(NULL),
77   fEnd(NULL),
78   fIndex(0x0)
79 {
80 // create an object to read digits from the given input file for the
81 // event with the given number
82
83   TDirectory* dir = gDirectory;
84   fFile = TFile::Open(fileName);
85   dir->cd();
86   if (!fFile || !fFile->IsOpen()) {
87     Error("AliRawReaderRoot", "could not open file %s", fileName);
88     fIsValid = kFALSE;
89     return;
90   }
91   TTree* tree = (TTree*) fFile->Get("RAW");
92   if (!tree) {
93     Error("AliRawReaderRoot", "no raw data tree found");
94     fIsValid = kFALSE;
95     return;
96   }
97   fBranch = tree->GetBranch("rawevent");
98   if (!fBranch) {
99     Error("AliRawReaderRoot", "no raw data branch found");
100     fIsValid = kFALSE;
101     return;
102   }
103
104   fBranch->SetAddress(&fEvent);
105   if (fEventIndex >= 0) {
106     if (fBranch->GetEntry(fEventIndex) <= 0) {
107       Error("AliRawReaderRoot", "no event with number %d found", fEventIndex);
108       fIsValid = kFALSE;
109       return;
110     }
111     fEventHeader = fEvent->GetHeader();
112   }
113 }
114
115 AliRawReaderRoot::AliRawReaderRoot(AliRawVEvent* event) :
116   fFile(NULL),
117   fBranch(NULL),
118   fEventIndex(-1),
119   fEvent(event),
120   fEventHeader(event->GetHeader()),
121   fSubEventIndex(0),
122   fSubEvent(NULL),
123   fEquipmentIndex(0),
124   fEquipment(NULL),
125   fRawData(NULL),
126   fPosition(NULL),
127   fEnd(NULL),
128   fIndex(0x0)
129 {
130 // create an object to read digits from the given raw event
131   if (!fEvent) fIsValid = kFALSE;
132 }
133
134 AliRawReaderRoot::AliRawReaderRoot(const AliRawReaderRoot& rawReader) :
135   AliRawReader(rawReader),
136   fFile(NULL),
137   fBranch(NULL),
138   fEventIndex(rawReader.fEventIndex),
139   fEvent(NULL),
140   fEventHeader(NULL),
141   fSubEventIndex(rawReader.fSubEventIndex),
142   fSubEvent(NULL),
143   fEquipmentIndex(rawReader.fEquipmentIndex),
144   fEquipment(NULL),
145   fRawData(NULL),
146   fPosition(NULL),
147   fEnd(NULL),
148   fIndex(0x0)
149 {
150 // copy constructor
151
152   if (rawReader.fFile) {
153     TDirectory* dir = gDirectory;
154     fFile = TFile::Open(rawReader.fFile->GetName());
155     dir->cd();
156     if (!fFile || !fFile->IsOpen()) {
157       Error("AliRawReaderRoot", "could not open file %s", 
158             rawReader.fFile->GetName());
159       fIsValid = kFALSE;
160       return;
161     }
162     TTree* tree = (TTree*) fFile->Get("RAW");
163     if (!tree) {
164       Error("AliRawReaderRoot", "no raw data tree found");
165       fIsValid = kFALSE;
166       return;
167     }
168     fBranch = tree->GetBranch("rawevent");
169     if (!fBranch) {
170       Error("AliRawReaderRoot", "no raw data branch found");
171       fIsValid = kFALSE;
172       return;
173     }
174
175     fBranch->SetAddress(&fEvent);
176     if (fEventIndex >= 0) {
177       if (fBranch->GetEntry(fEventIndex) <= 0) {
178         Error("AliRawReaderRoot", "no event with number %d found", 
179               fEventIndex);
180         fIsValid = kFALSE;
181         return;
182       }
183       fEventHeader = fEvent->GetHeader();
184     }
185   } else {
186     fEvent = rawReader.fEvent;
187     fEventHeader = rawReader.fEventHeader;
188   }
189
190   if (fSubEventIndex > 0) {
191     fSubEvent = fEvent->GetSubEvent(fSubEventIndex-1);
192     fEquipment = fSubEvent->GetEquipment(fEquipmentIndex);
193     fRawData = fEquipment->GetRawData();
194     fCount = 0;
195     fHeader = (AliRawDataHeader*) ((UChar_t*) fRawData->GetBuffer() + 
196       ((UChar_t*) rawReader.fHeader - 
197        (UChar_t*) rawReader.fRawData->GetBuffer()));
198     // Now check the version of the header
199     UChar_t version = 2;
200     if (fHeader) version = fHeader->GetVersion();
201     if (version==3) {
202       fHeader = NULL;
203       fHeaderV3 = (AliRawDataHeaderV3*) ((UChar_t*) fRawData->GetBuffer() + 
204       ((UChar_t*) rawReader.fHeaderV3 - 
205        (UChar_t*) rawReader.fRawData->GetBuffer()));
206     }
207     fPosition = (UChar_t*) fRawData->GetBuffer() + 
208       (rawReader.fPosition - (UChar_t*) rawReader.fRawData->GetBuffer());
209     fEnd = ((UChar_t*) fRawData->GetBuffer()) + fRawData->GetSize();
210   }
211 }
212
213 AliRawReaderRoot& AliRawReaderRoot::operator = (const AliRawReaderRoot& 
214                                                 rawReader)
215 {
216 // assignment operator
217
218   this->~AliRawReaderRoot();
219   new(this) AliRawReaderRoot(rawReader);
220   return *this;
221 }
222
223 AliRawReaderRoot::~AliRawReaderRoot()
224 {
225 // delete objects and close root file
226
227   if (fFile) {
228     if (fEvent) delete fEvent;
229     fFile->Close();
230     delete fFile;
231   }
232 }
233
234 const AliRawEventHeaderBase* AliRawReaderRoot::GetEventHeader() const
235 {
236   // Get the even header
237   // Return NULL in case of failure
238   return fEventHeader;
239 }
240
241 UInt_t AliRawReaderRoot::GetType() const
242 {
243 // get the type from the event header
244
245   if (!fEventHeader) return 0;
246   return fEventHeader->Get("Type");
247 }
248
249 UInt_t AliRawReaderRoot::GetRunNumber() const
250 {
251 // get the run number from the event header
252
253   if (!fEventHeader) return 0;
254   return fEventHeader->Get("RunNb");
255 }
256
257 const UInt_t* AliRawReaderRoot::GetEventId() const
258 {
259 // get the event id from the event header
260
261   if (!fEventHeader) return NULL;
262   return fEventHeader->GetP("Id");
263 }
264
265 const UInt_t* AliRawReaderRoot::GetTriggerPattern() const
266 {
267 // get the trigger pattern from the event header
268
269   if (!fEventHeader) return NULL;
270   return fEventHeader->GetP("TriggerPattern");
271 }
272
273 const UInt_t* AliRawReaderRoot::GetDetectorPattern() const
274 {
275 // get the detector pattern from the event header
276
277   if (!fEventHeader) return NULL;
278   return fEventHeader->GetP("DetectorPattern");
279 }
280
281 const UInt_t* AliRawReaderRoot::GetAttributes() const
282 {
283 // get the type attributes from the event header
284
285   if (!fEventHeader) return NULL;
286   return fEventHeader->GetP("TypeAttribute");
287 }
288
289 const UInt_t* AliRawReaderRoot::GetSubEventAttributes() const
290 {
291 // get the type attributes from the sub event header
292
293   if (!fSubEvent) return NULL;
294   return fSubEvent->GetHeader()->GetP("TypeAttribute");
295 }
296
297 UInt_t AliRawReaderRoot::GetLDCId() const
298 {
299 // get the LDC Id from the event header
300
301   if (!fEvent || !fSubEvent) return 0;
302   return fSubEvent->GetHeader()->Get("LdcId");
303 }
304
305 UInt_t AliRawReaderRoot::GetGDCId() const
306 {
307 // get the GDC Id from the event header
308
309   if (!fEventHeader) return 0;
310   return fEventHeader->Get("GdcId");
311 }
312
313 UInt_t AliRawReaderRoot::GetTimestamp() const
314 {
315   if (!fEventHeader) return 0;
316   return fEventHeader->Get("Timestamp");
317 }
318
319 Int_t AliRawReaderRoot::GetEquipmentSize() const
320 {
321 // get the size of the equipment
322
323   if (!fEvent || !fEquipment || !fEquipment->GetEquipmentHeader()) return 0;
324   return fEquipment->GetEquipmentHeader()->GetEquipmentSize();
325 }
326
327 Int_t AliRawReaderRoot::GetEquipmentType() const
328 {
329 // get the type from the equipment header
330
331   if (!fEvent || !fEquipment || !fEquipment->GetEquipmentHeader()) return -1;
332   return fEquipment->GetEquipmentHeader()->GetEquipmentType();
333 }
334
335 Int_t AliRawReaderRoot::GetEquipmentId() const
336 {
337 // get the ID from the equipment header
338
339   if (!fEvent || !fEquipment || !fEquipment->GetEquipmentHeader()) return -1;
340   return fEquipment->GetEquipmentHeader()->GetId();
341 }
342
343 const UInt_t* AliRawReaderRoot::GetEquipmentAttributes() const
344 {
345 // get the attributes from the equipment header
346
347   if (!fEvent || !fEquipment || !fEquipment->GetEquipmentHeader()) return NULL;
348   return fEquipment->GetEquipmentHeader()->GetTypeAttribute();
349 }
350
351 Int_t AliRawReaderRoot::GetEquipmentElementSize() const
352 {
353 // get the basic element size from the equipment header
354
355   if (!fEvent || !fEquipment || !fEquipment->GetEquipmentHeader()) return 0;
356   return fEquipment->GetEquipmentHeader()->GetBasicSizeType();
357 }
358
359 Int_t AliRawReaderRoot::GetEquipmentHeaderSize() const
360 {
361 // get the size of the equipment header (28 bytes by default)
362
363   if (!fEvent || !fEquipment || !fEquipment->GetEquipmentHeader()) return 0;
364   return fEquipment->GetEquipmentHeader()->HeaderSize();
365 }
366
367 // _________________________________________________________________________
368 void AliRawReaderRoot::SwapData(const void* inbuf, const void* outbuf, UInt_t size) {
369   // The method swaps the contents of the
370   // raw-data event header
371   UInt_t  intCount = (size+3)/sizeof(UInt_t);
372
373   UInt_t* buf = (UInt_t*) inbuf;    // temporary integers buffer
374   for (UInt_t i=0; i<intCount; i++, buf++) {
375       UInt_t value = SwapWord(*buf);
376       if (i==(intCount-1))
377          memcpy((UInt_t*)outbuf+i, &value, size%sizeof(UInt_t));
378       else
379          memcpy((UInt_t*)outbuf+i, &value, sizeof(UInt_t));
380   }
381 }
382 // _________________________________________________________________________
383
384 Bool_t AliRawReaderRoot::ReadHeader()
385 {
386 // read a data header at the current position
387 // returns kFALSE if the data header could not be read
388
389   fErrorCode = 0;
390   if (!fEvent) return kFALSE;
391
392   do {
393     // skip payload (if event was not selected)
394     if (fCount > 0) fPosition += fCount;
395
396     // get the first or the next equipment if at the end of an equipment
397     if (!fEquipment || (fPosition >= fEnd)) {
398
399       // get the first or the next sub event if at the end of a sub event
400       if (!fSubEvent || (fEquipmentIndex >= fSubEvent->GetNEquipments())) {
401
402         // check for end of event data
403         if (fSubEventIndex >= fEvent->GetNSubEvents()) return kFALSE;
404         fSubEvent = fEvent->GetSubEvent(fSubEventIndex++);
405
406         // check the magic word of the sub event
407         if (!fSubEvent->GetHeader()->IsValid()) {
408           Error("ReadHeader", "wrong magic number in sub event!");
409           fSubEvent->GetHeader()->Dump();
410           fErrorCode = kErrMagic;
411           return kFALSE;
412         }
413
414         fEquipmentIndex = 0;
415         fEquipment = NULL;
416         fRawData = NULL;
417       }
418
419       // get the next equipment and raw data
420       fCount = 0;
421       if (fEquipmentIndex >= fSubEvent->GetNEquipments()) {
422         fEquipment = NULL;
423         continue;
424       }
425       fEquipment = fSubEvent->GetEquipment(fEquipmentIndex++);
426       if (!fEquipment) continue;
427       if (!IsSelected()) {
428         fPosition = fEnd;
429         continue;
430       }
431       fRawData = fEquipment->GetRawData();
432       if (!fRawData) {
433         fPosition = fEnd;
434         continue;
435       }
436       fPosition = (UChar_t*) fRawData->GetBuffer();
437       fEnd = ((UChar_t*) fRawData->GetBuffer()) + fRawData->GetSize();
438     }
439
440     // continue with the next equipment if no data left in the payload
441     if (fPosition >= fEnd) continue;
442
443     if (fRequireHeader) {
444       // check that there are enough bytes left for the data header
445       if (fPosition + sizeof(AliRawDataHeader) > fEnd) {
446         Error("ReadHeader", "could not read data header!");
447         Warning("ReadHeader", "skipping %ld bytes", fEnd - fPosition);
448         fEquipment->GetEquipmentHeader()->Dump();
449         fCount = 0;
450         fPosition = fEnd;
451         fErrorCode = kErrNoDataHeader;
452         continue;
453       }
454
455       // "read" the data header
456       fHeader = (AliRawDataHeader*) fPosition;
457       // Now check the version of the header
458       UChar_t version = 2;
459       if (fHeader) version=fHeader->GetVersion();
460       if (version==2) {
461  #ifndef R__BYTESWAP
462         SwapData((void*) fHeader, (void*) fHeaderSwapped, sizeof(AliRawDataHeader));
463         fHeader=fHeaderSwapped;
464 #endif
465         if ((fPosition + fHeader->fSize) != fEnd) {
466           if (fHeader->fSize != 0xFFFFFFFF)
467             Warning("ReadHeader",
468                     "Equipment %d : raw data size found in the header is wrong (%d != %ld)! Using the equipment size instead !",
469                     fEquipment->GetEquipmentHeader()->GetId(),fHeader->fSize, fEnd - fPosition);
470           fHeader->fSize = fEnd - fPosition;
471         }
472         fPosition += sizeof(AliRawDataHeader);
473         fHeaderV3 = 0;
474       } else if (version==3) {
475         fHeaderV3 = (AliRawDataHeaderV3*) fPosition;
476 #ifndef R__BYTESWAP
477         SwapData((void*) fHeaderV3, (void*) fHeaderSwapped, sizeof(AliRawDataHeaderV3));
478         fHeaderV3=fHeaderSwappedV3;
479 #endif
480         if ((fPosition + fHeaderV3->fSize) != fEnd) {
481           if (fHeaderV3->fSize != 0xFFFFFFFF)
482             Warning("ReadHeader",
483                     "Equipment %d : raw data size found in the header is wrong (%d != %ld)! Using the equipment size instead !",
484                     fEquipment->GetEquipmentHeader()->GetId(),fHeader->fSize, fEnd - fPosition);
485           fHeaderV3->fSize = fEnd - fPosition;
486         }
487         fPosition += sizeof(AliRawDataHeaderV3);
488         fHeader = 0;
489       }
490     }
491
492     if (fHeader && (fHeader->fSize != 0xFFFFFFFF)) {
493       fCount = fHeader->fSize - sizeof(AliRawDataHeader);
494
495       // check consistency of data size in the header and in the sub event
496       if (fPosition + fCount > fEnd) {  
497         Error("ReadHeader", "size in data header exceeds event size!");
498         Warning("ReadHeader", "skipping %ld bytes", fEnd - fPosition);
499         fEquipment->GetEquipmentHeader()->Dump();
500         fCount = 0;
501         fPosition = fEnd;
502         fErrorCode = kErrSize;
503         continue;
504       }
505     } else if (fHeaderV3 && (fHeaderV3->fSize != 0xFFFFFFFF)) {
506       fCount = fHeaderV3->fSize - sizeof(AliRawDataHeaderV3);
507
508       // check consistency of data size in the header and in the sub event
509       if (fPosition + fCount > fEnd) {  
510         Error("ReadHeader", "size in data header exceeds event size!");
511         Warning("ReadHeader", "skipping %ld bytes", fEnd - fPosition);
512         fEquipment->GetEquipmentHeader()->Dump();
513         fCount = 0;
514         fPosition = fEnd;
515         fErrorCode = kErrSize;
516         continue;
517       }
518
519     } else {
520       fCount = fEnd - fPosition;
521     }
522
523   } while (!IsSelected());
524
525   return kTRUE;
526 }
527
528 Bool_t AliRawReaderRoot::ReadNextData(UChar_t*& data)
529 {
530 // reads the next payload at the current position
531 // returns kFALSE if the data could not be read
532
533   fErrorCode = 0;
534   while (fCount == 0) {
535     if (!ReadHeader()) return kFALSE;
536   }
537   data = fPosition;
538   fPosition += fCount;  
539   fCount = 0;
540   return kTRUE;
541 }
542
543 Bool_t AliRawReaderRoot::ReadNext(UChar_t* data, Int_t size)
544 {
545 // reads the next block of data at the current position
546 // returns kFALSE if the data could not be read
547
548   fErrorCode = 0;
549   if (fPosition + size > fEnd) {
550     Error("ReadNext", "could not read data!");
551     fErrorCode = kErrOutOfBounds;
552     return kFALSE;
553   }
554   memcpy(data, fPosition, size);
555
556   fPosition += size;
557   fCount -= size;
558   return kTRUE;
559 }
560
561
562 Bool_t AliRawReaderRoot::Reset()
563 {
564 // reset the current position to the beginning of the event
565
566   fSubEventIndex = 0;
567   fSubEvent = NULL;
568   fEquipmentIndex = 0;
569   fEquipment = NULL;
570   fRawData = NULL;
571   fHeader = NULL;
572   fHeaderV3 = NULL;
573
574   fCount = 0;
575   fPosition = fEnd = NULL;
576   return kTRUE;
577 }
578
579
580 Bool_t AliRawReaderRoot::NextEvent()
581 {
582 // go to the next event in the root file
583
584   if (!fBranch) return kFALSE;
585
586   // check if it uses order or not
587   if (fgUseOrder && !fIndex) MakeIndex();
588
589   do {
590     delete fEvent;
591     fEvent = NULL;
592     fEventHeader = NULL;
593     fBranch->SetAddress(&fEvent);
594     Int_t entryToGet = fEventIndex + 1;
595     if (fgUseOrder 
596         && fIndex 
597         && entryToGet<fBranch->GetEntries()
598         && entryToGet>-1 ) entryToGet = fIndex[entryToGet];
599     if (fBranch->GetEntry(entryToGet) <= 0)
600       return kFALSE;
601     fEventHeader = fEvent->GetHeader();
602     fEventIndex++;
603   } while (!IsEventSelected());
604   fEventNumber++;
605   return Reset();
606 }
607
608 Bool_t AliRawReaderRoot::RewindEvents()
609 {
610 // go back to the beginning of the root file
611
612   if (!fBranch) return kFALSE;
613
614   fEventIndex = -1;
615   delete fEvent;
616   fEvent = NULL;
617   fEventHeader = NULL;
618   fBranch->SetAddress(&fEvent);
619   fEventNumber = -1;
620   return Reset();
621 }
622
623 Bool_t  AliRawReaderRoot::GotoEvent(Int_t event)
624 {
625   // go to a particular event
626   // Uses the absolute event index inside the
627   // raw-data file
628
629   if (!fBranch) return kFALSE;
630
631   // check if it uses order or not
632   if (fgUseOrder && !fIndex) MakeIndex();
633
634   delete fEvent;
635   fEvent = NULL;
636   fEventHeader = NULL;
637   fBranch->SetAddress(&fEvent);
638   Int_t entryToGet = event;
639   if (fgUseOrder 
640       && fIndex 
641       && entryToGet<fBranch->GetEntries()
642       && entryToGet>-1 ) entryToGet = fIndex[entryToGet];
643   if (fBranch->GetEntry(entryToGet) <= 0)
644     return kFALSE;
645   fEventHeader = fEvent->GetHeader();
646   fEventIndex = event;
647   fEventNumber++;
648   return Reset();
649 }
650
651 Int_t AliRawReaderRoot::GetNumberOfEvents() const
652 {
653   // Get the total number of events in
654   // the raw-data tree
655
656   if (!fBranch) return -1;
657
658   return fBranch->GetEntries();
659 }
660
661 Int_t AliRawReaderRoot::CheckData() const
662 {
663 // check the consistency of the data
664
665   if (!fEvent) return 0;
666
667   AliRawVEvent* subEvent = NULL;
668   Int_t subEventIndex = 0;
669   AliRawVEquipment* equipment = NULL;
670   Int_t equipmentIndex = 0;
671   UChar_t* position = 0;
672   UChar_t* end = 0;
673   Int_t result = 0;
674
675   while (kTRUE) {
676     // get the first or the next sub event if at the end of an equipment
677     if (!subEvent || (equipmentIndex >= subEvent->GetNEquipments())) {
678
679       // check for end of event data
680       if (subEventIndex >= fEvent->GetNSubEvents()) return result;
681       subEvent = fEvent->GetSubEvent(subEventIndex++);
682
683       // check the magic word of the sub event
684       if (!fSubEvent->GetHeader()->IsValid()) {
685         result |= kErrMagic;
686         return result;
687       }
688
689       equipmentIndex = 0;
690     }
691
692     // get the next equipment and raw data
693     if (equipmentIndex >= subEvent->GetNEquipments()) {
694       equipment = NULL;
695       continue;
696     }
697     equipment = subEvent->GetEquipment(equipmentIndex++);
698     if (!equipment) continue;
699     AliRawData* rawData = equipment->GetRawData();
700     if (!rawData) continue;
701     position = (UChar_t*) rawData->GetBuffer();
702     end = ((UChar_t*) rawData->GetBuffer()) + rawData->GetSize();
703
704     // continue with the next sub event if no data left in the payload
705     if (position >= end) continue;
706
707     if (fRequireHeader) {
708     // check that there are enough bytes left for the data header
709       if (position + sizeof(AliRawDataHeader) > end) {
710         result |= kErrNoDataHeader;
711         continue;
712       }
713
714       // check consistency of data size in the header and in the equipment
715       AliRawDataHeader* header = (AliRawDataHeader*) position;
716       if ((position + header->fSize) != end) {
717         if (header->fSize != 0xFFFFFFFF)
718           Warning("ReadHeader",
719                   "Equipment %d : raw data size found in the header is wrong (%d != %ld)! Using the equipment size instead !",
720                   equipment->GetEquipmentHeader()->GetId(),header->fSize, end - position);
721         header->fSize = end - position;
722         result |= kErrSize;
723       }
724     }
725     position = end;
726   };
727
728   return result;
729 }
730
731 AliRawReader* AliRawReaderRoot::CloneSingleEvent() const
732 {
733   // Clones the current event and
734   // creates raw-reader for the cloned event
735   // Can be used in order to make asynchronious
736   // access to the current raw data within
737   // several threads (online event display/reco)
738
739   if (fEvent) {
740     // Root formatted raw data
741     AliRawVEvent *gdcRootEvent = (AliRawVEvent*)fEvent->Clone();
742     for (Int_t ldcCounter=0; ldcCounter < gdcRootEvent->GetNSubEvents(); ldcCounter++) {
743       AliRawVEvent *ldcRootEvent = gdcRootEvent->GetSubEvent(ldcCounter);
744       AliRawVEvent *subEvent = fEvent->GetSubEvent(ldcCounter);
745       for (Int_t eqCounter=0; eqCounter < ldcRootEvent->GetNEquipments(); eqCounter++) {
746         AliRawVEquipment *equipment=ldcRootEvent->GetEquipment(eqCounter);
747         AliRawVEquipment *eq = subEvent->GetEquipment(eqCounter);
748         equipment->CloneRawData(eq->GetRawData());
749       }
750     }
751     // Reset original event and newly
752     // produced one
753     gdcRootEvent->GetSubEvent(-1);
754     fEvent->GetSubEvent(-1);
755     return new AliRawReaderRoot(gdcRootEvent);
756   }
757   return NULL;
758 }
759
760 void AliRawReaderRoot::MakeIndex() {
761   // Make index
762   if (fBranch) {
763     TTree * rawTree = fBranch->GetTree();
764     if (rawTree) {
765       rawTree->BuildIndex("-fEvtHdrs[0].fSize"); // Minus sign to get largest first
766       TTreeIndex * treeInd = (TTreeIndex*)rawTree->GetTreeIndex();
767       if (treeInd) fIndex = treeInd->GetIndex();
768     }
769   }
770 }