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