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