]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDrawData.cxx
Missing DATE event types are added to the base raw data header
[u/mrichter/AliRoot.git] / TRD / AliTRDrawData.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 //  TRD raw data conversion class                                            //
21 //                                                                           //
22 ///////////////////////////////////////////////////////////////////////////////
23
24 #include <Riostream.h>
25
26 #include "AliTRDrawData.h"
27 #include "AliTRDdigitsManager.h"
28 #include "AliTRDgeometry.h"
29 #include "AliTRDdataArrayI.h"
30 #include "AliTRDRawStream.h"
31 #include "AliRawDataHeader.h"
32 #include "AliRawReader.h"
33 #include "AliTRDCommonParam.h"
34 #include "AliTRDcalibDB.h"
35 #include "AliDAQ.h"
36
37 ClassImp(AliTRDrawData)
38
39 //_____________________________________________________________________________
40 AliTRDrawData::AliTRDrawData():TObject()
41 {
42   //
43   // Default constructor
44   //
45
46   fDebug         = 0;
47
48 }
49
50 //_____________________________________________________________________________
51 AliTRDrawData::AliTRDrawData(const AliTRDrawData &r):TObject()
52 {
53   //
54   // AliTRDrawData copy constructor
55   //
56
57   ((AliTRDrawData &) r).Copy(*this);
58
59 }
60
61 //_____________________________________________________________________________
62 AliTRDrawData::~AliTRDrawData()
63 {
64   //
65   // Destructor
66   //
67
68 }
69
70 //_____________________________________________________________________________
71 AliTRDrawData &AliTRDrawData::operator=(const AliTRDrawData &r)
72 {
73   //
74   // Assignment operator
75   //
76
77   if (this != &r) ((AliTRDrawData &) r).Copy(*this);
78   return *this;
79
80 }
81
82 //_____________________________________________________________________________
83 void AliTRDrawData::Copy(TObject &r) const
84 {
85   //
86   // Copy function
87   //
88
89   ((AliTRDrawData &) r).fDebug         = fDebug;
90
91 }
92
93 //_____________________________________________________________________________
94 Bool_t AliTRDrawData::Digits2Raw(TTree *digitsTree)
95 {
96   //
97   // Convert the digits to raw data byte stream. The output is written
98   // into the the binary files TRD_<DDL number>.ddl.
99   //
100   // The pseudo raw data format is currently defined like this:
101   //
102   //          DDL data header
103   //
104   //          Subevent (= single chamber) header (8 bytes)
105   //                  FLAG
106   //                  Detector number (2 bytes)
107   //                  Number of data bytes (2 bytes)
108   //                  Number of pads with data (2 bytes)
109   //                  1 empty byte
110   //
111   //          Data bank
112   //
113
114   const Int_t kNumberOfDDLs         = AliDAQ::NumberOfDdls("TRD");
115   const Int_t kSubeventHeaderLength = 8;
116   const Int_t kSubeventDummyFlag    = 0xBB;
117   Int_t       headerSubevent[3];
118
119   ofstream     **outputFile = new ofstream* [kNumberOfDDLs];
120   UInt_t        *bHPosition = new UInt_t    [kNumberOfDDLs];
121   Int_t         *ntotalbyte = new Int_t     [kNumberOfDDLs];
122   Int_t          nbyte = 0;
123   Int_t          npads = 0;
124   unsigned char *bytePtr;
125   unsigned char *headerPtr;
126
127   AliTRDdigitsManager* digitsManager = new AliTRDdigitsManager();
128   digitsManager->SetDebug(fDebug);
129
130   // Read in the digit arrays
131   if (!digitsManager->ReadDigits(digitsTree)) {
132     delete digitsManager;
133     return kFALSE;
134   }
135
136   AliTRDgeometry   *geo = new AliTRDgeometry();
137   AliTRDdataArrayI *digits;
138
139   AliTRDCommonParam* commonParam = AliTRDCommonParam::Instance();
140   if (!commonParam)
141   {
142     printf("<AliTRDrawData::Digits2Raw> ");
143     printf("Could not get common params\n");
144     return 0;
145   }
146   
147   AliTRDcalibDB* calibration = AliTRDcalibDB::Instance();
148   if (!calibration)
149   {
150     printf("<AliTRDdigitizer::Digits2Raw> ");
151     printf("Could not get calibration object\n");
152     return kFALSE;
153   }
154     
155   // the event header
156   AliRawDataHeader header;
157
158   // Open the output files
159   for (Int_t iDDL = 0; iDDL < kNumberOfDDLs; iDDL++) {
160     char name[20];
161     strcpy(name,AliDAQ::DdlFileName("TRD",iDDL));
162 #ifndef __DECCXX
163     outputFile[iDDL] = new ofstream(name, ios::binary);
164 #else
165     outputFile[iDDL] = new ofstream(name);
166 #endif
167
168     // Write a dummy data header
169     bHPosition[iDDL] = outputFile[iDDL]->tellp();
170     outputFile[iDDL]->write((char*)(&header),sizeof(header));
171     ntotalbyte[iDDL] = 0;
172   }
173
174   // Loop through all detectors
175   for (Int_t det = 0; det < AliTRDgeometry::Ndet(); det++) {
176
177     Int_t cham      = geo->GetChamber(det);
178     Int_t plan      = geo->GetPlane(det);
179     Int_t sect      = geo->GetSector(det);
180     Int_t rowMax    = commonParam->GetRowMax(plan,cham,sect);
181     Int_t colMax    = commonParam->GetColMax(plan);
182     Int_t timeTotal = calibration->GetNumberOfTimeBins();
183     Int_t bufferMax = rowMax*colMax*timeTotal;
184     Int_t *buffer   = new Int_t[bufferMax];
185
186     npads   = 0;
187     nbyte   = 0;
188     bytePtr = (unsigned char *) buffer;
189
190     Int_t iDDL = sect;
191
192     // Get the digits array
193     digits = digitsManager->GetDigits(det);
194     digits->Expand();
195
196     // Loop through the detector pixel
197     for (Int_t col = 0; col < colMax; col++) {
198       for (Int_t row = 0; row < rowMax; row++) {
199
200         // Check whether data exists for this pad
201         Bool_t dataflag = kFALSE;
202         for (Int_t time = 0; time < timeTotal; time++) {
203           Int_t data = digits->GetDataUnchecked(row,col,time);
204           if (data) {
205             dataflag = kTRUE;
206             break;
207           }
208         }
209
210         if (dataflag) {
211
212           npads++;
213
214           // The pad row number
215           *bytePtr++ = row + 1;
216           // The pad column number
217           *bytePtr++ = col + 1;
218           nbyte += 2;
219
220           Int_t nzero = 0;
221           for (Int_t time = 0; time < timeTotal; time++) {
222
223             Int_t data = digits->GetDataUnchecked(row,col,time);
224
225             if (!data) {
226               nzero++;
227               if ((nzero ==       256) || 
228                   (time  == timeTotal-1)) {
229                 *bytePtr++ = 0;
230                 *bytePtr++ = nzero-1;
231                 nbyte += 2;
232                 nzero  = 0;
233               }
234             }
235             else {
236               if (nzero) {
237                 *bytePtr++ = 0;
238                 *bytePtr++ = nzero-1;
239                 nbyte += 2;
240                 nzero  = 0;
241               }
242               // High byte (MSB always set)
243               *bytePtr++ = ((data >> 8) | 128);
244               // Low byte
245               *bytePtr++ = (data & 0xff);
246               nbyte += 2;
247             }
248
249           }
250
251         }
252
253       }
254
255     }
256
257     // Fill the end of the buffer with zeros
258     while (nbyte % 4) {  
259       *bytePtr++ = 0;
260       nbyte++;
261     }
262
263     if (fDebug > 1) {
264       Info("Digits2Raw","det = %d, nbyte = %d (%d)",det,nbyte,bufferMax);
265     }
266
267     // Write the subevent header
268     bytePtr    = (unsigned char *) headerSubevent;
269     headerPtr  = bytePtr;
270     *bytePtr++ = kSubeventDummyFlag;
271     *bytePtr++ = (det   & 0xff);
272     *bytePtr++ = (det   >> 8);
273     *bytePtr++ = (nbyte & 0xff);
274     *bytePtr++ = (nbyte >> 8);
275     *bytePtr++ = (nbyte >> 16);
276     *bytePtr++ = (npads & 0xff);
277     *bytePtr++ = (npads >> 8);
278     outputFile[iDDL]->write((char*)headerPtr,kSubeventHeaderLength);
279
280     // Write the buffer to the file
281     bytePtr = (unsigned char *) buffer;
282     outputFile[iDDL]->write((char*)bytePtr,nbyte);
283
284     ntotalbyte[iDDL] += nbyte + kSubeventHeaderLength;
285
286     delete buffer;
287
288   }
289
290   if (fDebug) {
291     for (Int_t iDDL = 0; iDDL < kNumberOfDDLs; iDDL++) {
292       Info("Digits2Raw","Total size: DDL %d = %d",iDDL,ntotalbyte[iDDL]);
293     }
294   }
295
296   // Update the data headers and close the output files
297   for (Int_t iDDL = 0; iDDL < kNumberOfDDLs; iDDL++) {
298     header.fSize = UInt_t(outputFile[iDDL]->tellp()) - bHPosition[iDDL];
299     header.SetAttribute(0);  // valid data
300     outputFile[iDDL]->seekp(bHPosition[iDDL]);
301     outputFile[iDDL]->write((char*)(&header),sizeof(header));
302
303     outputFile[iDDL]->close();
304     delete outputFile[iDDL];
305   }
306
307   delete geo;
308   delete digitsManager;
309
310   delete [] outputFile;
311   delete [] bHPosition;
312   delete [] ntotalbyte;
313
314
315
316
317   return kTRUE;
318
319 }
320
321 //_____________________________________________________________________________
322 AliTRDdigitsManager* AliTRDrawData::Raw2Digits(AliRawReader* rawReader)
323 {
324   //
325   // Read the raw data digits and put them into the returned digits manager
326   //
327
328   AliTRDdataArrayI *digits    = 0;
329   AliTRDdataArrayI *track0    = 0;
330   AliTRDdataArrayI *track1    = 0;
331   AliTRDdataArrayI *track2    = 0; 
332
333   AliTRDgeometry *geo = new AliTRDgeometry();
334
335   AliTRDCommonParam* commonParam = AliTRDCommonParam::Instance();
336   if (!commonParam)
337   {
338     printf("<AliTRDrawData::Raw2Digits> ");
339     printf("Could not get common params\n");
340     return 0;
341   }
342     
343   AliTRDcalibDB* calibration = AliTRDcalibDB::Instance();
344   if (!calibration)
345   {
346     printf("<AliTRDdigitizer::Raw2Digits> ");
347     printf("Could not get calibration object\n");
348     return 0;
349   }
350
351   // Create the digits manager
352   AliTRDdigitsManager* digitsManager = new AliTRDdigitsManager();
353   digitsManager->SetDebug(fDebug);
354   digitsManager->CreateArrays();
355
356   AliTRDRawStream input(rawReader);
357
358   // Loop through the digits
359   while (input.Next()) {
360
361     Int_t det    = input.GetDetector();
362     Int_t npads  = input.GetNPads();
363
364     if (input.IsNewDetector()) {
365
366       if (digits) digits->Compress(1,0);
367       if (track0) track0->Compress(1,0);
368       if (track1) track1->Compress(1,0);
369       if (track2) track2->Compress(1,0);
370
371       if (fDebug > 2) {
372         Info("Raw2Digits","Subevent header:");
373         Info("Raw2Digits","\tdet   = %d",det);
374         Info("Raw2Digits","\tnpads = %d",npads);
375       }      
376
377       // Create the data buffer
378       Int_t cham      = geo->GetChamber(det);
379       Int_t plan      = geo->GetPlane(det);
380       Int_t sect      = geo->GetSector(det);
381       Int_t rowMax    = commonParam->GetRowMax(plan,cham,sect);
382       Int_t colMax    = commonParam->GetColMax(plan);
383       Int_t timeTotal = calibration->GetNumberOfTimeBins();
384
385       // Add a container for the digits of this detector
386       digits = digitsManager->GetDigits(det);
387       track0 = digitsManager->GetDictionary(det,0);
388       track1 = digitsManager->GetDictionary(det,1);
389       track2 = digitsManager->GetDictionary(det,2);
390       // Allocate memory space for the digits buffer
391       if (digits->GetNtime() == 0) {
392         digits->Allocate(rowMax,colMax,timeTotal);
393         track0->Allocate(rowMax,colMax,timeTotal);
394         track1->Allocate(rowMax,colMax,timeTotal);
395         track2->Allocate(rowMax,colMax,timeTotal);
396       }
397
398     } 
399
400     digits->SetDataUnchecked(input.GetRow(),input.GetColumn(),
401                              input.GetTime(),input.GetSignal());
402     track0->SetDataUnchecked(input.GetRow(),input.GetColumn(),
403                              input.GetTime(),               -1);
404     track1->SetDataUnchecked(input.GetRow(),input.GetColumn(),
405                              input.GetTime(),               -1);
406     track2->SetDataUnchecked(input.GetRow(),input.GetColumn(),
407                              input.GetTime(),               -1);
408   }
409
410   if (digits) digits->Compress(1,0);
411   if (track0) track0->Compress(1,0);
412   if (track1) track1->Compress(1,0);
413   if (track2) track2->Compress(1,0);
414
415   delete geo;
416
417   return digitsManager;
418
419 }