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