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