]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PMD/AliPMDDDLRawData.cxx
change in mapping
[u/mrichter/AliRoot.git] / PMD / AliPMDDDLRawData.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 #include <Riostream.h>
17 #include <TClonesArray.h>
18 #include <TTree.h>
19 #include <TBranch.h>
20 #include <TMath.h>
21
22 #include "AliRawDataHeader.h"
23 #include "AliPMDdigit.h"
24 #include "AliPMDRawStream.h"
25 #include "AliPMDDDLRawData.h"
26
27
28 ClassImp(AliPMDDDLRawData)
29
30 AliPMDDDLRawData::AliPMDDDLRawData():
31   fDigits(new TClonesArray("AliPMDdigit", 1000))
32 {
33   // Default Constructor
34   //
35
36 }
37 //____________________________________________________________________________
38
39 AliPMDDDLRawData::~AliPMDDDLRawData()
40 {
41   // Default Destructor
42   //
43
44 }
45
46 //____________________________________________________________________________
47 void AliPMDDDLRawData::WritePMDRawData(TTree *treeD)
48 {
49   // write digits into raw data format
50
51   ofstream outfile;
52
53   TBranch *branch = treeD->GetBranch("PMDDigit");
54   if (!branch) return;
55   branch->SetAddress(&fDigits);  
56   
57   //  Int_t   nmodules = (Int_t) treeD->GetEntries();
58   //  cout << " nmodules = " << nmodules << endl;
59
60   const Int_t kSize         = 4608;
61   const Int_t kDDL          = 6;
62   Int_t modulePerDDL        = 0;
63
64
65   AliRawDataHeader header;
66   UInt_t sizeRawData = 0;
67   Int_t  totword     = 0;
68
69   UInt_t buffer[kSize];
70
71   Char_t filename[80];
72
73   Int_t mmodule = 0;
74   for(Int_t iddl = 0; iddl < kDDL; iddl++)
75     {
76       sprintf(filename,"PMD_%d.ddl",AliPMDRawStream::kDDLOffset+iddl);
77 #ifndef __DECCXX
78       outfile.open(filename,ios::binary);
79 #else
80       outfile.open(filename);
81 #endif
82       
83       if (iddl < 4)
84         {
85           modulePerDDL = 6;
86           mmodule = 6*iddl;
87         }
88       else if (iddl == 4)
89         {
90           modulePerDDL = 12;
91           mmodule = 24;
92         }
93       else if (iddl == 5)
94         {
95           modulePerDDL = 12;
96           mmodule = 30;
97         }
98
99
100       // Write the Dummy Data Header into the file
101       Int_t bHPosition = outfile.tellp();
102       outfile.write((char*)(&header),sizeof(header));
103
104
105       for(Int_t ium = 0; ium < modulePerDDL; ium++)
106         {
107           
108           if (iddl == 4 && ium == 6) mmodule = 36;
109           if (iddl == 5 && ium == 6) mmodule = 42;
110
111           for (Int_t i = 0; i < kSize; i++)
112             {
113               buffer[i] = 0;
114             }
115           // Extract energy deposition per cell and pack it
116           // in a 32-bit word and returns all the total words
117           // per one unit-module
118           
119           GetUMDigitsData(treeD, mmodule, ium, iddl, totword, buffer);
120           
121           outfile.write((char*)buffer,totword*sizeof(UInt_t));
122
123           mmodule++;
124
125         }
126
127       // Write real data header
128       // take the pointer to the beginning of the data header
129       // write the total number of words per ddl and bring the
130       // pointer to the current file position and close it
131       UInt_t cFPosition = outfile.tellp();
132       sizeRawData = cFPosition - bHPosition - sizeof(header);
133       header.fSize = cFPosition - bHPosition;
134       header.SetAttribute(0);  // valid data
135       outfile.seekp(bHPosition);
136       outfile.write((char*)(&header),sizeof(header));
137       outfile.seekp(cFPosition);
138
139       outfile.close();
140     } // DDL Loop over
141
142
143 }
144 //____________________________________________________________________________
145 void AliPMDDDLRawData::GetUMDigitsData(TTree *treeD, Int_t imodule, Int_t ium,
146                                        Int_t ddlno, Int_t & totword,
147                                        UInt_t *buffer)
148 {
149
150   UInt_t dataword, baseword;
151
152   UInt_t mcmno, chno;
153   UInt_t adc;
154   Int_t  irownew = 0;
155   Int_t  icolnew = 0;
156   Int_t  irownew1 = 0;
157   Int_t  icolnew1 = 0;
158   Int_t  det, smn, irow, icol;
159
160
161   treeD->GetEntry(imodule); 
162   Int_t nentries = fDigits->GetLast();
163   totword = nentries+1;
164
165   for (Int_t ient = 0; ient < totword; ient++)
166     {
167       fPMDdigit = (AliPMDdigit*)fDigits->UncheckedAt(ient);
168       
169       det    = fPMDdigit->GetDetector();
170       smn    = fPMDdigit->GetSMNumber();
171       irow   = fPMDdigit->GetRow();
172       icol   = fPMDdigit->GetColumn();
173       adc    = (UInt_t) fPMDdigit->GetADC();
174
175       //      cout << " imodule = " << imodule << " smn = " << smn << endl;
176
177       if(smn < 12)
178         {
179           irownew1 = icol;
180           icolnew1 = irow;
181         }
182       else if( smn >= 12 && smn < 24)
183         {
184           irownew1 = irow;
185           icolnew1 = icol;
186         }
187
188       if(smn < 6)
189         {
190           irownew = 95 - irownew1;
191           icolnew = icolnew1;
192         }
193       else if(smn >= 6 && smn < 12)
194         {
195           irownew = irownew1;
196           icolnew = 47 - icolnew1;
197         }
198       else if(smn >= 12 && smn < 18)
199         {
200           irownew = 47 - irownew1;
201           icolnew = icolnew1;
202         }
203       else if(smn >= 18 && smn < 24)
204         {
205           irownew = irownew1;
206           icolnew = 95 - icolnew1;
207         }
208       
209
210
211       GetMCMCh(ddlno, ium, irownew, icolnew, mcmno, chno);
212
213       baseword = 0;
214       dataword = adc;
215       PackWord(0, 11,dataword,baseword);
216       dataword = chno;
217       PackWord(12,17,dataword,baseword);
218       dataword = mcmno;
219       PackWord(18,28,dataword,baseword);
220       dataword = 0;
221       PackWord(29,29,dataword,baseword);
222       dataword = 0;
223       PackWord(30,30,dataword,baseword);
224       dataword = 1;
225       PackWord(31,31,dataword,baseword);
226       
227       buffer[ient] = baseword;
228       
229     }
230
231
232 }
233 //____________________________________________________________________________
234 void AliPMDDDLRawData::GetMCMCh(Int_t ddlno, Int_t um,
235                                 Int_t row, Int_t col,
236                                 UInt_t &mcmno, UInt_t &chno)
237 {
238   // This part will be modified once the final track layout on the PCB is
239   // designed. This will only change the coordinate of the individual cell
240
241   static const UInt_t kCh[16][4] = { {12, 13, 18, 19},
242                                      {11, 15, 17, 21},
243                                      {10, 14, 16, 22},
244                                      {8, 9, 20, 23},
245                                      {7, 4, 25, 24},
246                                      {6, 0, 30, 26},
247                                      {5, 1, 31, 27},
248                                      {3, 2, 29, 28},
249                                      {44, 45, 50, 51},
250                                      {43, 47, 49, 53},
251                                      {42, 46, 48, 54},
252                                      {40, 41, 52, 55},
253                                      {39, 36, 57, 56},
254                                      {38, 32, 62, 58},
255                                      {37, 33, 63, 59},
256                                      {35, 34, 61, 60} };
257   
258   if (ddlno == 0 || ddlno == 1)
259     {
260       // PRE plane, SU Mod = 1, 2
261       Int_t irownew = row%16;
262       Int_t icolnew = col%4;
263       Int_t irowdiv = row/16;
264       Int_t icoldiv = col/4;
265
266       mcmno = 72*um + 12*irowdiv + icoldiv;
267       chno  = kCh[irownew][icolnew];
268     }
269   else if (ddlno == 2 || ddlno == 3)
270     {
271       // PRE plane,  SU Mod = 3, 4
272       Int_t irownew = row%16;
273       Int_t icolnew = col%4;
274       Int_t irowdiv = row/16;
275       Int_t icoldiv = col/4;
276       
277       mcmno = 72*um + 24*irowdiv + icoldiv;
278       chno  = kCh[irownew][icolnew];
279     }
280   else if (ddlno == 4 || ddlno == 5)
281     {
282       // CPV plane,  SU Mod = 1, 3 : ddl = 4
283       // CPV plane,  SU Mod = 2, 4 : ddl = 5
284
285       Int_t irownew = row%16;
286       Int_t icolnew = col%4;
287       Int_t irowdiv = row/16;
288       Int_t icoldiv = col/4;
289
290       if(um < 6)
291         {
292           mcmno = 72*um + 12*irowdiv + icoldiv;
293         }
294       else if(um >= 6)
295         {
296           mcmno = 72*um + 24*irowdiv + icoldiv;
297         }
298       chno  = kCh[irownew][icolnew];
299     }
300
301 }
302 //____________________________________________________________________________
303
304 void AliPMDDDLRawData::PackWord(UInt_t startbit, UInt_t stopbit, 
305                                 UInt_t dataword, UInt_t &packedword)
306 {
307   UInt_t bitLength  = stopbit - startbit + 1;
308   UInt_t bitContent = (UInt_t) (TMath::Power(2,bitLength) - 1);
309   if(bitContent < dataword)
310     {
311       cout << " *** ERROR *** bitContent is less than the dataword" << endl;
312       return;
313     }
314   UInt_t packedBits = 0;
315   if (packedword != 0)
316     packedBits = (UInt_t) (TMath::Log(packedword)/TMath::Log(2));
317
318   UInt_t counter;
319   if (packedBits <= stopbit)
320     {
321       counter   = 31 - stopbit;
322     }
323   else
324     {
325       counter   = 31 - packedBits;
326     }
327   UInt_t dummyword = 0xFFFFFFFF;
328   dummyword >>= counter;
329   UInt_t lword = dataword << startbit;
330   UInt_t nword = lword | packedword;
331   packedword = dummyword & nword;
332
333
334 }
335 //____________________________________________________________________________
336 void AliPMDDDLRawData::UnpackWord(UInt_t startbit, UInt_t stopbit, 
337                                 UInt_t &dataword, UInt_t packedword)
338 {
339   UInt_t bitLength  = stopbit - startbit + 1;
340   UInt_t bitContent = (UInt_t) (TMath::Power(2,bitLength) - 1);
341   bitContent <<= startbit;
342   dataword = packedword & bitContent;
343   dataword >>= startbit;
344
345 }
346 //____________________________________________________________________________
347