]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RAW/AliRawReaderRoot.cxx
Minor modification related to the new CASTOR stager
[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 /// AliRawEvent.
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 "AliRawReaderRoot.h"
34 #include "AliRawEvent.h"
35 #include "AliRawEventHeader.h"
36 #include "AliRawEquipment.h"
37 #include "AliRawEquipmentHeader.h"
38 #include "AliRawData.h"
39
40
41 ClassImp(AliRawReaderRoot)
42
43
44 AliRawReaderRoot::AliRawReaderRoot(const char* fileName, Int_t eventNumber) :
45   fFile(NULL),
46   fBranch(NULL),
47   fEventIndex(eventNumber),
48   fEvent(NULL),
49   fSubEventIndex(0),
50   fSubEvent(NULL),
51   fEquipmentIndex(0),
52   fEquipment(NULL),
53   fRawData(NULL),
54   fPosition(NULL),
55   fEnd(NULL)
56 {
57 // create an object to read digits from the given input file for the
58 // event with the given number
59
60   TDirectory* dir = gDirectory;
61   fFile = TFile::Open(fileName);
62   dir->cd();
63   if (!fFile || !fFile->IsOpen()) {
64     Error("AliRawReaderRoot", "could not open file %s", fileName);
65     return;
66   }
67   TTree* tree = (TTree*) fFile->Get("RAW");
68   if (!tree) {
69     Error("AliRawReaderRoot", "no raw data tree found");
70     return;
71   }
72   fBranch = tree->GetBranch("rawevent");
73   if (!fBranch) {
74     Error("AliRawReaderRoot", "no raw data branch found");
75     return;
76   }
77
78   fEvent = new AliRawEvent;
79   fBranch->SetAddress(&fEvent);
80   if (fEventIndex >= 0) {
81     if (fBranch->GetEntry(fEventIndex) <= 0) {
82       Error("AliRawReaderRoot", "no event with number %d found", fEventIndex);
83       return;
84     }
85   }
86 }
87
88 AliRawReaderRoot::AliRawReaderRoot(AliRawEvent* event) :
89   fFile(NULL),
90   fBranch(NULL),
91   fEventIndex(-1),
92   fEvent(event),
93   fSubEventIndex(0),
94   fSubEvent(NULL),
95   fEquipmentIndex(0),
96   fEquipment(NULL),
97   fRawData(NULL),
98   fPosition(NULL),
99   fEnd(NULL)
100 {
101 // create an object to read digits from the given raw event
102
103 }
104
105 AliRawReaderRoot::AliRawReaderRoot(const AliRawReaderRoot& rawReader) :
106   AliRawReader(rawReader),
107   fFile(NULL),
108   fBranch(NULL),
109   fEventIndex(rawReader.fEventIndex),
110   fEvent(NULL),
111   fSubEventIndex(rawReader.fSubEventIndex),
112   fSubEvent(NULL),
113   fEquipmentIndex(rawReader.fEquipmentIndex),
114   fEquipment(NULL),
115   fRawData(NULL),
116   fPosition(NULL),
117   fEnd(NULL)
118 {
119 // copy constructor
120
121   if (rawReader.fFile) {
122     TDirectory* dir = gDirectory;
123     fFile = TFile::Open(rawReader.fFile->GetName());
124     dir->cd();
125     if (!fFile || !fFile->IsOpen()) {
126       Error("AliRawReaderRoot", "could not open file %s", 
127             rawReader.fFile->GetName());
128       return;
129     }
130     TTree* tree = (TTree*) fFile->Get("RAW");
131     if (!tree) {
132       Error("AliRawReaderRoot", "no raw data tree found");
133       return;
134     }
135     fBranch = tree->GetBranch("rawevent");
136     if (!fBranch) {
137       Error("AliRawReaderRoot", "no raw data branch found");
138       return;
139     }
140
141     fEvent = new AliRawEvent;
142     fBranch->SetAddress(&fEvent);
143     if (fEventIndex >= 0) {
144       if (fBranch->GetEntry(fEventIndex) <= 0) {
145         Error("AliRawReaderRoot", "no event with number %d found", 
146               fEventIndex);
147         return;
148       }
149     }
150   } else {
151     fEvent = rawReader.fEvent;
152   }
153
154   if (fSubEventIndex > 0) {
155     fSubEvent = fEvent->GetSubEvent(fSubEventIndex-1);
156     fEquipment = fSubEvent->GetEquipment(fEquipmentIndex);
157     fRawData = fEquipment->GetRawData();
158       fCount = 0;
159     fHeader = (AliRawDataHeader*) ((UChar_t*) fRawData->GetBuffer() + 
160       ((UChar_t*) rawReader.fHeader - 
161        (UChar_t*) rawReader.fRawData->GetBuffer()));
162     fPosition = (UChar_t*) fRawData->GetBuffer() + 
163       (rawReader.fPosition - (UChar_t*) rawReader.fRawData->GetBuffer());
164     fEnd = ((UChar_t*) fRawData->GetBuffer()) + fRawData->GetSize();
165   }
166 }
167
168 AliRawReaderRoot& AliRawReaderRoot::operator = (const AliRawReaderRoot& 
169                                                 rawReader)
170 {
171 // assignment operator
172
173   this->~AliRawReaderRoot();
174   new(this) AliRawReaderRoot(rawReader);
175   return *this;
176 }
177
178 AliRawReaderRoot::~AliRawReaderRoot()
179 {
180 // delete objects and close root file
181
182   if (fFile) {
183     if (fEvent) delete fEvent;
184     fFile->Close();
185     delete fFile;
186   }
187 }
188
189
190 UInt_t AliRawReaderRoot::GetType() const
191 {
192 // get the type from the event header
193
194   if (!fEvent) return 0;
195   return fEvent->GetHeader()->GetType();
196 }
197
198 UInt_t AliRawReaderRoot::GetRunNumber() const
199 {
200 // get the run number from the event header
201
202   if (!fEvent) return 0;
203   return fEvent->GetHeader()->GetRunNumber();
204 }
205
206 const UInt_t* AliRawReaderRoot::GetEventId() const
207 {
208 // get the event id from the event header
209
210   if (!fEvent) return NULL;
211   return fEvent->GetHeader()->GetId();
212 }
213
214 const UInt_t* AliRawReaderRoot::GetTriggerPattern() const
215 {
216 // get the trigger pattern from the event header
217
218   if (!fEvent) return NULL;
219   return fEvent->GetHeader()->GetTriggerPattern();
220 }
221
222 const UInt_t* AliRawReaderRoot::GetDetectorPattern() const
223 {
224 // get the detector pattern from the event header
225
226   if (!fEvent) return NULL;
227   return fEvent->GetHeader()->GetDetectorPattern();
228 }
229
230 const UInt_t* AliRawReaderRoot::GetAttributes() const
231 {
232 // get the type attributes from the event header
233
234   if (!fEvent) return NULL;
235   return fEvent->GetHeader()->GetTypeAttribute();
236 }
237
238 const UInt_t* AliRawReaderRoot::GetSubEventAttributes() const
239 {
240 // get the type attributes from the sub event header
241
242   if (!fSubEvent) return NULL;
243   return fSubEvent->GetHeader()->GetTypeAttribute();
244 }
245
246 UInt_t AliRawReaderRoot::GetLDCId() const
247 {
248 // get the LDC Id from the event header
249
250   if (!fEvent || !fSubEvent) return 0;
251   return fSubEvent->GetHeader()->GetLDCId();
252 }
253
254 UInt_t AliRawReaderRoot::GetGDCId() const
255 {
256 // get the GDC Id from the event header
257
258   if (!fEvent) return 0;
259   return fEvent->GetHeader()->GetGDCId();
260 }
261
262
263 Int_t AliRawReaderRoot::GetEquipmentSize() const
264 {
265 // get the size of the equipment
266
267   if (!fEvent || !fEquipment || !fEquipment->GetEquipmentHeader()) return 0;
268   return fEquipment->GetEquipmentHeader()->GetEquipmentSize();
269 }
270
271 Int_t AliRawReaderRoot::GetEquipmentType() const
272 {
273 // get the type from the equipment header
274
275   if (!fEvent || !fEquipment || !fEquipment->GetEquipmentHeader()) return -1;
276   return fEquipment->GetEquipmentHeader()->GetEquipmentType();
277 }
278
279 Int_t AliRawReaderRoot::GetEquipmentId() const
280 {
281 // get the ID from the equipment header
282
283   if (!fEvent || !fEquipment || !fEquipment->GetEquipmentHeader()) return -1;
284   return fEquipment->GetEquipmentHeader()->GetId();
285 }
286
287 const UInt_t* AliRawReaderRoot::GetEquipmentAttributes() const
288 {
289 // get the attributes from the equipment header
290
291   if (!fEvent || !fEquipment || !fEquipment->GetEquipmentHeader()) return NULL;
292   return fEquipment->GetEquipmentHeader()->GetTypeAttribute();
293 }
294
295 Int_t AliRawReaderRoot::GetEquipmentElementSize() const
296 {
297 // get the basic element size from the equipment header
298
299   if (!fEvent || !fEquipment || !fEquipment->GetEquipmentHeader()) return 0;
300   return fEquipment->GetEquipmentHeader()->GetBasicSizeType();
301 }
302
303
304 Bool_t AliRawReaderRoot::ReadHeader()
305 {
306 // read a data header at the current position
307 // returns kFALSE if the data header could not be read
308
309   fErrorCode = 0;
310   if (!fEvent) return kFALSE;
311
312   do {
313     // skip payload (if event was not selected)
314     if (fCount > 0) fPosition += fCount;
315
316     // get the first or the next equipment if at the end of an equipment
317     if (!fEquipment || (fPosition >= fEnd)) {
318
319       // get the first or the next sub event if at the end of a sub event
320       if (!fSubEvent || (fEquipmentIndex >= fSubEvent->GetNEquipments())) {
321
322         // check for end of event data
323         if (fSubEventIndex >= fEvent->GetNSubEvents()) return kFALSE;
324         fSubEvent = fEvent->GetSubEvent(fSubEventIndex++);
325
326         // check the magic word of the sub event
327         if (!fSubEvent->GetHeader()->IsValid()) {
328           Error("ReadHeader", "wrong magic number in sub event!");
329           fSubEvent->GetHeader()->Dump();
330           fErrorCode = kErrMagic;
331           return kFALSE;
332         }
333
334         fEquipmentIndex = 0;
335         fEquipment = NULL;
336         fRawData = NULL;
337       }
338
339       // get the next equipment and raw data
340       fCount = 0;
341       fEquipment = fSubEvent->GetEquipment(fEquipmentIndex++);
342       if (!fEquipment) continue;
343       fRawData = fEquipment->GetRawData();
344       if (!fRawData) {
345         fPosition = fEnd;
346         continue;
347       }
348       fPosition = (UChar_t*) fRawData->GetBuffer();
349       fEnd = ((UChar_t*) fRawData->GetBuffer()) + fRawData->GetSize();
350     }
351
352     // continue with the next equipment if no data left in the payload
353     if (fPosition >= fEnd) continue;
354
355     // check that there are enough bytes left for the data header
356     if (fPosition + sizeof(AliRawDataHeader) > fEnd) {
357       Error("ReadHeader", "could not read data header!");
358       Warning("ReadHeader", "skipping %d bytes", fEnd - fPosition);
359       fEquipment->GetEquipmentHeader()->Dump();
360       fCount = 0;
361       fPosition = fEnd;
362       fErrorCode = kErrNoDataHeader;
363       continue;
364     }
365
366     // "read" the data header
367     fHeader = (AliRawDataHeader*) fPosition;
368     fPosition += sizeof(AliRawDataHeader);
369     if (fHeader->fSize != 0xFFFFFFFF) {
370       fCount = fHeader->fSize - sizeof(AliRawDataHeader);
371     } else {
372       fCount = fEnd - fPosition;
373     }
374
375     // check consistency of data size in the header and in the sub event
376     if (fPosition + fCount > fEnd) {  
377       Error("ReadHeader", "size in data header exceeds event size!");
378       Warning("ReadHeader", "skipping %d bytes", fEnd - fPosition);
379       fEquipment->GetEquipmentHeader()->Dump();
380       fCount = 0;
381       fPosition = fEnd;
382       fErrorCode = kErrSize;
383       continue;
384     }
385
386   } while (!IsSelected());
387
388   return kTRUE;
389 }
390
391 Bool_t AliRawReaderRoot::ReadNextData(UChar_t*& data)
392 {
393 // reads the next payload at the current position
394 // returns kFALSE if the data could not be read
395
396   fErrorCode = 0;
397   while (fCount == 0) {
398     if (!ReadHeader()) return kFALSE;
399   }
400   data = fPosition;
401   fPosition += fCount;  
402   fCount = 0;
403   return kTRUE;
404 }
405
406 Bool_t AliRawReaderRoot::ReadNext(UChar_t* data, Int_t size)
407 {
408 // reads the next block of data at the current position
409 // returns kFALSE if the data could not be read
410
411   fErrorCode = 0;
412   if (fPosition + size > fEnd) {
413     Error("ReadNext", "could not read data!");
414     fErrorCode = kErrOutOfBounds;
415     return kFALSE;
416   }
417   memcpy(data, fPosition, size);
418   fPosition += size;
419   fCount -= size;
420   return kTRUE;
421 }
422
423
424 Bool_t AliRawReaderRoot::Reset()
425 {
426 // reset the current position to the beginning of the event
427
428   fSubEventIndex = 0;
429   fSubEvent = NULL;
430   fEquipmentIndex = 0;
431   fEquipment = NULL;
432   fRawData = NULL;
433   fHeader = NULL;
434
435   fCount = 0;
436   fPosition = fEnd = NULL;
437   return kTRUE;
438 }
439
440
441 Bool_t AliRawReaderRoot::NextEvent()
442 {
443 // go to the next event in the root file
444
445   if (!fFile) return kFALSE;
446
447   do {
448     if (fBranch->GetEntry(fEventIndex+1) <= 0) {
449       delete fEvent;
450       fEvent = new AliRawEvent;
451       fBranch->SetAddress(&fEvent);
452       return kFALSE;
453     }
454     fEventIndex++;
455   } while (!IsEventSelected());
456   return Reset();
457 }
458
459 Bool_t AliRawReaderRoot::RewindEvents()
460 {
461 // go back to the beginning of the root file
462
463   if (!fFile) return kFALSE;
464
465   fEventIndex = -1;
466   delete fEvent;
467   fEvent = new AliRawEvent;
468   fBranch->SetAddress(&fEvent);
469   return Reset();
470 }
471
472
473 Int_t AliRawReaderRoot::CheckData() const
474 {
475 // check the consistency of the data
476
477   if (!fEvent) return 0;
478
479   AliRawEvent* subEvent = NULL;
480   Int_t subEventIndex = 0;
481   AliRawEquipment* equipment = NULL;
482   Int_t equipmentIndex = 0;
483   UChar_t* position = 0;
484   UChar_t* end = 0;
485   Int_t result = 0;
486
487   while (kTRUE) {
488     // get the first or the next sub event if at the end of an equipment
489     if (!subEvent || (equipmentIndex >= subEvent->GetNEquipments())) {
490
491       // check for end of event data
492       if (subEventIndex >= fEvent->GetNSubEvents()) return result;
493       subEvent = fEvent->GetSubEvent(subEventIndex++);
494
495       // check the magic word of the sub event
496       if (!fSubEvent->GetHeader()->IsValid()) {
497         result |= kErrMagic;
498         return result;
499       }
500
501       equipmentIndex = 0;
502     }
503
504     // get the next equipment and raw data
505     equipment = subEvent->GetEquipment(equipmentIndex++);
506     if (!equipment) continue;
507     AliRawData* rawData = equipment->GetRawData();
508     if (!rawData) continue;
509     position = (UChar_t*) rawData->GetBuffer();
510     end = ((UChar_t*) rawData->GetBuffer()) + rawData->GetSize();
511
512     // continue with the next sub event if no data left in the payload
513     if (position >= end) continue;
514
515     // check that there are enough bytes left for the data header
516     if (position + sizeof(AliRawDataHeader) > end) {
517       result |= kErrNoDataHeader;
518       continue;
519     }
520
521     // check consistency of data size in the header and in the equipment
522     AliRawDataHeader* header = (AliRawDataHeader*) position;
523     if (header->fSize != 0xFFFFFFFF) {
524       if (position + header->fSize > end) {
525         result |= kErrSize;
526       }
527     }
528   };
529
530   return result;
531 }