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