]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PMD/AliPMDDDLRawData.cxx
PMD DDL raw data
[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 "AliPMDdigit.h"
23 #include "AliPMDDDLRawData.h"
24
25
26 ClassImp(AliPMDDDLRawData)
27
28 AliPMDDDLRawData::AliPMDDDLRawData():
29   fDigits(new TClonesArray("AliPMDdigit", 1000))
30 {
31   // Default Constructor
32   //
33 }
34 //____________________________________________________________________________
35
36 AliPMDDDLRawData::~AliPMDDDLRawData()
37 {
38   // Default Destructor
39   //
40
41 }
42
43 //____________________________________________________________________________
44 void AliPMDDDLRawData::WritePMDRawData(TTree *treeD, Int_t evtno)
45 {
46   // write digits into raw data format
47
48   ofstream outfile;
49
50   TBranch *branch = treeD->GetBranch("PMDDigit");
51   if (!branch) return;
52   branch->SetAddress(&fDigits);  
53   
54   Int_t   nmodules = (Int_t) treeD->GetEntries();
55   //  cout << " nmodules = " << nmodules << endl;
56
57   const Int_t kSize         = 4608;
58   const Int_t kDDL          = 6;
59   Int_t modulePerDDL        = 0;
60
61
62   UInt_t sizeRawData = 0;
63   UInt_t magicWord   = 0x123456;
64   UInt_t detectorID  = 0;
65   UInt_t ddlID       = 0;
66   Int_t  flag        = 0;
67   Int_t  version     = 1;
68   UInt_t mHSize      = 3*sizeof(UInt_t);
69   Int_t  totword     = 0;
70
71   UInt_t miniHeader[3];
72   UInt_t buffer[kSize];
73
74   Char_t filename[80];
75
76   Int_t mmodule = 0;
77   for(Int_t iddl = 0; iddl < kDDL; iddl++)
78     {
79       sprintf(filename,"Ev%dPMDddl%d.dat",evtno,iddl);
80 #ifndef __DECCXX
81       outfile.open(filename,ios::binary);
82 #else
83       outfile.open(filename);
84 #endif
85       
86       if (iddl < 4)
87         {
88           modulePerDDL = 6;
89           mmodule = 6*iddl;
90           detectorID = 0;
91         }
92       else if (iddl == 4)
93         {
94           modulePerDDL = 12;
95           mmodule = 24;
96           detectorID = 1;
97         }
98       else if (iddl == 5)
99         {
100           modulePerDDL = 12;
101           mmodule = 36;
102           detectorID = 1;
103         }
104
105       miniHeader[0] = sizeRawData;
106       PackWord(0,23,magicWord,miniHeader[1]);
107       PackWord(24,31,detectorID,miniHeader[1]);
108       ddlID = iddl;
109       PackWord(0,15,ddlID,miniHeader[2]);
110       PackWord(16,23,flag,miniHeader[2]);
111       PackWord(24,31,version,miniHeader[2]);
112
113
114       // Write the Dummy Mini Header into the file
115       Int_t bHPosition = outfile.tellp();
116       outfile.write((char*)(miniHeader),mHSize);
117
118
119       for(Int_t ium = 0; ium < modulePerDDL; ium++)
120         {
121           
122           for (Int_t i = 0; i < kSize; i++)
123             {
124               buffer[i] = 0;
125             }
126           // Extract energy deposition per cell and pack it
127           // in a 32-bit word and returns all the total words
128           // per one unit-module
129           
130           GetUMDigitsData(treeD, mmodule, ium, iddl, totword, buffer);
131           
132           outfile.write((char*)buffer,totword*sizeof(UInt_t));
133
134           mmodule++;
135
136         }
137
138       // Write real mini header
139       // take the pointer to the beginning of the mini header
140       // write the total number of words per ddl and bring the
141       // pointer to the current file position and close it
142       UInt_t cFPosition = outfile.tellp();
143       sizeRawData = cFPosition - bHPosition - mHSize;
144       outfile.seekp(bHPosition);
145       outfile.write((char*)(&sizeRawData),sizeof(UInt_t));
146       outfile.seekp(cFPosition);
147
148       outfile.close();
149     } // DDL Loop over
150
151
152 }
153 //____________________________________________________________________________
154
155 void AliPMDDDLRawData::ReadPMDRawData(Int_t evtno)
156 {
157   // reads the raw data
158   ifstream infile;
159
160   Char_t filename[80];
161
162   const Int_t kDDL = 6;
163   Int_t imodule = 0;
164
165   for(Int_t iddl = 0; iddl < kDDL; iddl++)
166     {
167       sprintf(filename,"Ev%dPMDddl%d.dat",evtno,iddl);
168 #ifndef __DECCXX
169       infile.open(filename,ios::binary);
170 #else
171       infile.open(filename);
172 #endif
173       
174       Int_t  ium;
175       Int_t  irow;
176       Int_t  icol;
177       UInt_t baseword;
178       UInt_t miniHeader[3];
179       UInt_t sizeRawData;
180       UInt_t mHSize = 3*sizeof(UInt_t);
181
182       infile.read((char*)(miniHeader), mHSize);
183       
184       sizeRawData = miniHeader[0];
185       Int_t totword = sizeRawData/4;
186       
187       UInt_t adc;
188       UInt_t chno;
189       UInt_t mcmno;
190       
191       for (Int_t iword = 0; iword < totword; iword++)
192         {
193           infile.read((char*)(&baseword),sizeof(UInt_t));
194           
195           UnpackWord(0,11,adc,baseword);
196           UnpackWord(12,17,chno,baseword);
197           UnpackWord(18,28,mcmno,baseword);
198           
199           GetRowCol(iddl, mcmno, chno, ium, irow, icol);
200           if (iddl < 4)
201             {
202               imodule = iddl*6 + ium;
203             }
204           else if (iddl == 4)
205             {
206               imodule = 24 + ium;
207             }
208           else if (iddl == 5)
209             {
210               imodule = 36 + ium;
211             }
212
213           fprintf(fpw3,"%d %d %d %d %d %d\n",imodule,irow,icol,mcmno,chno,adc);
214           
215         }
216
217     } // ddl loop
218 }
219 //____________________________________________________________________________
220 void AliPMDDDLRawData::GetUMDigitsData(TTree *treeD, Int_t imodule, Int_t ium,
221                                        Int_t ddlno, Int_t & totword,
222                                        UInt_t *buffer)
223 {
224
225   UInt_t dataword, baseword;
226
227   UInt_t mcmno, chno;
228   UInt_t adc;
229   Int_t  irownew = 0;
230   Int_t  icolnew = 0;
231   Int_t  det, smn, irow, icol;
232
233
234   treeD->GetEntry(imodule); 
235   Int_t nentries = fDigits->GetLast();
236   totword = nentries+1;
237
238   for (Int_t ient = 0; ient < totword; ient++)
239     {
240       fPMDdigit = (AliPMDdigit*)fDigits->UncheckedAt(ient);
241       
242       det    = fPMDdigit->GetDetector();
243       smn    = fPMDdigit->GetSMNumber();
244       irow   = fPMDdigit->GetRow();
245       icol   = fPMDdigit->GetColumn();
246       adc    = (UInt_t) fPMDdigit->GetADC();
247
248       if(smn < 12)
249         {
250           irownew = icol;
251           icolnew = irow;
252         }
253       else if( smn >= 12 && smn < 24)
254         {
255           irownew = irow;
256           icolnew = icol;
257         }
258
259
260
261       GetMCMCh(ddlno, ium, irownew, icolnew, mcmno, chno);
262
263       fprintf(fpw2,"%d %d %d %d %d %d\n",imodule,irownew,icolnew,mcmno,chno,adc);
264
265       baseword = 0;
266       dataword = adc;
267       PackWord(0, 11,dataword,baseword);
268       dataword = chno;
269       PackWord(12,17,dataword,baseword);
270       dataword = mcmno;
271       PackWord(18,28,dataword,baseword);
272       dataword = 0;
273       PackWord(29,29,dataword,baseword);
274       dataword = 0;
275       PackWord(30,30,dataword,baseword);
276       dataword = 1;
277       PackWord(31,31,dataword,baseword);
278       
279       buffer[ient] = baseword;
280       
281
282
283       //      fprintf(fpw1,"%d %d %d %d %d\n",det,smn,irow,icol,adc);
284       
285     }
286
287
288 }
289 //____________________________________________________________________________
290 void AliPMDDDLRawData::GetMCMCh(Int_t ddlno, Int_t um,
291                                 Int_t row, Int_t col,
292                                 UInt_t &mcmno, UInt_t &chno)
293 {
294   // This part will be modified once the final track layout on the PCB is
295   // designed. This will only change the coordinate of the individual cell
296   UInt_t ch[16][4] = { {3, 2, 29, 28},
297                        {5, 1, 31, 27},
298                        {6, 0, 30, 26},
299                        {7, 4, 25, 24},
300                        {8, 9, 20, 23},
301                        {10, 14, 16, 22},
302                        {11, 15, 17, 21},
303                        {12, 13, 18, 19},
304                        {35, 34, 61, 60},
305                        {37, 33, 63, 59},
306                        {38, 32, 62, 58},
307                        {39, 36, 57, 56},
308                        {40, 41, 52, 55},
309                        {42, 46, 48, 54},
310                        {43, 47, 49, 53},
311                        {44, 45, 50, 51} };
312   
313   if (ddlno == 0 || ddlno == 1)
314     {
315       // PRE plane, SU Mod = 1, 2
316       Int_t irownew = row%16;
317       Int_t icolnew = col%4;
318       Int_t irowdiv = row/16;
319       Int_t icoldiv = col/4;
320
321       mcmno = 72*um + 6*icoldiv + irowdiv;
322       chno  = ch[irownew][icolnew];
323     }
324   else if (ddlno == 2 || ddlno == 3)
325     {
326       // PRE plane,  SU Mod = 3, 4
327       Int_t irownew = row%16;
328       Int_t icolnew = col%4;
329       Int_t irowdiv = row/16;
330       Int_t icoldiv = col/4;
331       
332       mcmno = 72*um + 3*icoldiv + irowdiv;
333       chno  = ch[irownew][icolnew];
334     }
335   else if (ddlno == 4)
336     {
337       // CPV plane,  SU Mod = 1, 2
338       Int_t irownew = row%16;
339       Int_t icolnew = col%4;
340       Int_t irowdiv = row/16;
341       Int_t icoldiv = col/4;
342
343       mcmno = 72*um + 6*icoldiv + irowdiv;
344       chno  = ch[irownew][icolnew];
345     }
346   else if (ddlno == 5)
347     {
348       // CPV plane,  SU Mod = 3, 4
349       Int_t irownew = row%16;
350       Int_t icolnew = col%4;
351       Int_t irowdiv = row/16;
352       Int_t icoldiv = col/4;
353       
354       mcmno = 72*um + 3*icoldiv + irowdiv;
355       chno  = ch[irownew][icolnew];
356     }
357
358 }
359 //____________________________________________________________________________
360
361 void AliPMDDDLRawData::GetRowCol(Int_t ddlno, UInt_t mcmno, UInt_t chno,
362                                  Int_t &um, Int_t &row, Int_t &col)
363 {
364   UInt_t ch[64] = { 9, 5, 1, 0, 13, 4, 8, 12,
365                     16, 17, 20, 24, 28, 29, 21, 25,
366                     22, 26, 30, 31, 18, 27, 23, 19,
367                     15, 14, 11, 7, 3, 2, 10, 6,
368                     41, 37, 33, 32, 45, 36, 40, 44,
369                     48, 49, 52, 56, 60, 61, 53, 57,
370                     54, 58, 62, 63, 50, 59, 55, 51,
371                     47, 46, 43, 39, 35, 34, 42, 38 };
372
373   if (ddlno == 1 || ddlno == 2)
374     {
375       um  = mcmno/72;
376       Int_t mcmnonew = mcmno - 72*um;
377       Int_t rowcol  = ch[chno];
378       Int_t irownew = rowcol/4;
379       Int_t icolnew = rowcol%4;
380       
381       Int_t remmcm  = mcmnonew%6;
382       Int_t divmcm  = mcmnonew/6;
383       
384       row = 16*remmcm + irownew;
385       col =  4*divmcm + icolnew;
386     }
387   else   if (ddlno == 3 || ddlno == 4)
388     {
389       um  = mcmno/72;
390       Int_t mcmnonew = mcmno - 72*um;
391       Int_t rowcol  = ch[chno];
392       Int_t irownew = rowcol/4;
393       Int_t icolnew = rowcol%4;
394       
395       Int_t remmcm  = mcmnonew%3;
396       Int_t divmcm  = mcmnonew/3;
397       
398       row = 16*remmcm + irownew;
399       col =  4*divmcm + icolnew;
400     }
401   else if (ddlno == 4)
402     {
403       um  = mcmno/144;
404       Int_t mcmnonew = mcmno - 72*um;
405       Int_t rowcol  = ch[chno];
406       Int_t irownew = rowcol/4;
407       Int_t icolnew = rowcol%4;
408       
409       Int_t remmcm  = mcmnonew%6;
410       Int_t divmcm  = mcmnonew/6;
411       
412       row = 16*remmcm + irownew;
413       col =  4*divmcm + icolnew;
414     }
415   else if (ddlno == 5)
416     {
417       um  = mcmno/144;
418       Int_t mcmnonew = mcmno - 72*um;
419       Int_t rowcol  = ch[chno];
420       Int_t irownew = rowcol/4;
421       Int_t icolnew = rowcol%4;
422       
423       Int_t remmcm  = mcmnonew%3;
424       Int_t divmcm  = mcmnonew/3;
425       
426       row = 16*remmcm + irownew;
427       col =  4*divmcm + icolnew;
428     }
429
430
431 }
432 //____________________________________________________________________________
433
434 void AliPMDDDLRawData::PackWord(UInt_t startbit, UInt_t stopbit, 
435                                 UInt_t dataword, UInt_t &packedword)
436 {
437   UInt_t bitLength  = stopbit - startbit + 1;
438   UInt_t bitContent = (UInt_t) (TMath::Power(2,bitLength) - 1);
439   if(bitContent < dataword)
440     {
441       cout << " *** ERROR *** bitContent is less than the dataword" << endl;
442       return;
443     }
444   UInt_t packedBits = 0;
445   if (packedword != 0)
446     packedBits = (UInt_t) (TMath::Log(packedword)/TMath::Log(2));
447
448   UInt_t counter;
449   if (packedBits <= stopbit)
450     {
451       counter   = 31 - stopbit;
452     }
453   else
454     {
455       counter   = 31 - packedBits;
456     }
457   UInt_t dummyword = 0xFFFFFFFF;
458   dummyword >>= counter;
459   UInt_t lword = dataword << startbit;
460   UInt_t nword = lword | packedword;
461   packedword = dummyword & nword;
462
463
464 }
465 //____________________________________________________________________________
466 void AliPMDDDLRawData::UnpackWord(UInt_t startbit, UInt_t stopbit, 
467                                 UInt_t &dataword, UInt_t packedword)
468 {
469   UInt_t bitLength  = stopbit - startbit + 1;
470   UInt_t bitContent = (UInt_t) (TMath::Power(2,bitLength) - 1);
471   bitContent <<= startbit;
472   dataword = packedword & bitContent;
473   dataword >>= startbit;
474
475 }
476 //____________________________________________________________________________
477