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