]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PMD/AliPMDRawStream.cxx
DP:Misalignment of CPV added
[u/mrichter/AliRoot.git] / PMD / AliPMDRawStream.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 /// This class provides access to PMD digits in raw data.
21 ///
22 /// It loops over all PMD digits in the raw data given by the AliRawReader.
23 /// The Next method goes to the next digit. If there are no digits left
24 /// it returns kFALSE.
25 /// Several getters provide information about the current digit.
26 ///
27 ///////////////////////////////////////////////////////////////////////////////
28
29 #include <Riostream.h>
30 #include <TObjArray.h>
31 #include <TString.h>
32 #include <TSystem.h>
33
34 #include "AliLog.h"
35 #include "AliPMDBlockHeader.h"
36 #include "AliPMDDspHeader.h"
37 #include "AliPMDPatchBusHeader.h"
38 #include "AliPMDddldata.h"
39 #include "AliPMDRawStream.h"
40 #include "AliRawReader.h"
41
42 ClassImp(AliPMDRawStream)
43
44
45 //_____________________________________________________________________________
46 AliPMDRawStream::AliPMDRawStream(AliRawReader* rawReader) :
47   fRawReader(rawReader)
48 {
49 // create an object to read PMD raw digits
50
51   fRawReader->Select("PMD");
52 }
53
54 //_____________________________________________________________________________
55 AliPMDRawStream::AliPMDRawStream(const AliPMDRawStream& stream) :
56   TObject(stream),
57   fRawReader(NULL)
58 {
59 // copy constructor
60
61   AliFatal("Copy constructor not implemented");
62 }
63
64 //_____________________________________________________________________________
65 AliPMDRawStream& AliPMDRawStream::operator = (const AliPMDRawStream& 
66                                               /* stream */)
67 {
68 // assignment operator
69
70   AliFatal("operator = assignment operator not implemented");
71   return *this;
72 }
73
74 //_____________________________________________________________________________
75 AliPMDRawStream::~AliPMDRawStream()
76 {
77 // destructor
78
79 }
80
81
82 //_____________________________________________________________________________
83
84 Bool_t AliPMDRawStream::DdlData(Int_t indexDDL, TObjArray *pmdddlcont)
85 {
86 // read the next raw digit
87 // returns kFALSE if there is no digit left
88
89   AliPMDddldata *pmdddldata;
90
91   if (!fRawReader->ReadHeader()) return kFALSE;
92   Int_t  iddl  = fRawReader->GetDDLID();
93   Int_t dataSize = fRawReader->GetDataSize();
94   Int_t totaldataword = dataSize/4;
95
96   if (dataSize <= 0) return kFALSE;
97   if (indexDDL != iddl)
98     {
99       AliError("Mismatch in the DDL index");
100       return kFALSE;
101     }
102
103   UInt_t *buffer;
104   buffer = new UInt_t[totaldataword];
105   UInt_t data;
106   for (Int_t i = 0; i < totaldataword; i++)
107     {
108       fRawReader->ReadNextInt(data);
109       buffer[i] = data;
110     }
111
112   // --- Open the mapping file
113
114   TString fileName(gSystem->Getenv("ALICE_ROOT"));
115   if(iddl == 0)
116     {
117       fileName += "/PMD/PMD_Mapping_ddl0.dat";
118     }
119   else if(iddl == 1)
120     {
121       fileName += "/PMD/PMD_Mapping_ddl1.dat";
122     }
123   else if(iddl == 2)
124     {
125       fileName += "/PMD/PMD_Mapping_ddl2.dat";
126     }
127   else if(iddl == 3)
128     {
129       fileName += "/PMD/PMD_Mapping_ddl3.dat";
130     }
131   else if(iddl == 4)
132     {
133       fileName += "/PMD/PMD_Mapping_ddl4.dat";
134     }
135   else if(iddl == 5)
136     {
137       fileName += "/PMD/PMD_Mapping_ddl5.dat";
138     }
139
140   ifstream infile;
141   infile.open(fileName.Data(), ios::in); // ascii file
142   if(!infile)
143     AliError(Form("Could not read the mapping file for DDL No = %d",iddl));
144   
145   Int_t modulePerDDL = 0;
146   if (iddl < 4)
147     {
148       modulePerDDL = 6;
149     }
150   else if (iddl == 4 || iddl == 5)
151     {
152       modulePerDDL = 12;
153     }
154
155   const Int_t kNPatchBus = 50;
156   
157   Int_t modno, totPatchBus, bPatchBus, ePatchBus;
158   Int_t ibus, totmcm, rows, rowe, cols, cole;
159   Int_t moduleNo[kNPatchBus], mcmperBus[kNPatchBus];
160   Int_t startRowBus[kNPatchBus], endRowBus[kNPatchBus];
161   Int_t startColBus[kNPatchBus], endColBus[kNPatchBus];
162
163
164   for (Int_t ibus = 0; ibus < kNPatchBus; ibus++)
165     {
166       mcmperBus[ibus]   = -1;
167       startRowBus[ibus] = -1;
168       endRowBus[ibus]   = -1;
169       startColBus[ibus] = -1;
170       endColBus[ibus]   = -1;
171     }
172
173
174   for (Int_t im = 0; im < modulePerDDL; im++)
175     {
176       infile >> modno;
177       infile >> totPatchBus >> bPatchBus >> ePatchBus;
178       
179       for(Int_t i=0; i<totPatchBus; i++)
180         {
181           infile >> ibus >> totmcm >> rows >> rowe >> cols >> cole;
182
183           moduleNo[ibus]     = modno;
184           mcmperBus[ibus]    = totmcm;
185           startRowBus[ibus]  = rows;
186           endRowBus[ibus]    = rowe;
187           startColBus[ibus]  = cols;
188           endColBus[ibus]    = cole;
189         }
190     }
191
192
193
194   infile.close();
195
196
197   AliPMDBlockHeader    blockHeader;
198   AliPMDDspHeader      dspHeader;
199   AliPMDPatchBusHeader pbusHeader;
200
201   const Int_t kblHLen   = blockHeader.GetHeaderLength();
202   const Int_t kdspHLen  = dspHeader.GetHeaderLength();
203   const Int_t kpbusHLen = pbusHeader.GetHeaderLength();
204   
205   Int_t parity;
206   Int_t idet, ismn;
207   Int_t irow = -1;
208   Int_t icol = -1;
209
210   Int_t blHeaderWord[8];
211   Int_t dspHeaderWord[10];
212   Int_t pbusHeaderWord[4];
213
214   Int_t ilowLimit = 0;
215   Int_t iuppLimit = 0;
216
217   Int_t blRawDataLength = 0;
218   Int_t iwordcount = 0;
219
220
221   for (Int_t iblock = 0; iblock < 2; iblock++)
222     {
223       ilowLimit = iuppLimit;
224       iuppLimit = ilowLimit + kblHLen;
225
226
227       for (Int_t i = ilowLimit; i < iuppLimit; i++)
228         {
229           blHeaderWord[i-ilowLimit] = (Int_t) buffer[i];
230         }
231
232       blockHeader.SetHeader(blHeaderWord);
233
234       blRawDataLength = blockHeader.GetRawDataLength();
235
236       for (Int_t idsp = 0; idsp < 5; idsp++)
237         {
238           ilowLimit = iuppLimit;
239           iuppLimit = ilowLimit + kdspHLen;
240
241           for (Int_t i = ilowLimit; i < iuppLimit; i++)
242             {
243               iwordcount++;
244               dspHeaderWord[i-ilowLimit] = (Int_t) buffer[i];
245             }
246           dspHeader.SetHeader(dspHeaderWord);
247
248           for (Int_t ibus = 0; ibus < 5; ibus++)
249             {
250
251               ilowLimit = iuppLimit;
252               iuppLimit = ilowLimit + kpbusHLen;
253
254               for (Int_t i = ilowLimit; i < iuppLimit; i++)
255                 {
256                   iwordcount++;
257                   pbusHeaderWord[i-ilowLimit] = (Int_t) buffer[i];
258                 }
259               pbusHeader.SetHeader(pbusHeaderWord);
260               Int_t rawdatalength = pbusHeader.GetRawDataLength();
261               Int_t pbusid = pbusHeader.GetPatchBusId();
262
263               ilowLimit = iuppLimit;
264               iuppLimit = ilowLimit + rawdatalength;
265
266               Int_t imodule = moduleNo[pbusid];
267
268
269               for (Int_t iword = ilowLimit; iword < iuppLimit; iword++)
270                 {
271                   iwordcount++;
272                   data = buffer[iword];
273
274                   Int_t isig =  data & 0x0FFF;
275                   Int_t ich  = (data >> 12) & 0x003F;
276                   Int_t imcm = (data >> 18) & 0x07FF;
277                   Int_t ibit = (data >> 31) & 0x0001;
278                   parity = ComputeParity(data);
279                   if (ibit != parity)
280                     {
281                       AliWarning("ComputeParity:: Parity Error");
282                     }
283                   GetRowCol(iddl, pbusid, imcm, ich, 
284                             startRowBus, endRowBus,
285                             startColBus, endColBus,
286                             irow, icol);
287
288                   ConvertDDL2SMN(iddl, imodule, ismn, idet);
289                   TransformH2S(ismn, irow, icol);
290
291                   pmdddldata = new AliPMDddldata();
292
293                   pmdddldata->SetDetector(idet);
294                   pmdddldata->SetSMN(ismn);
295                   pmdddldata->SetModule(imodule);
296                   pmdddldata->SetPatchBusId(pbusid);
297                   pmdddldata->SetMCM(imcm);
298                   pmdddldata->SetChannel(ich);
299                   pmdddldata->SetRow(irow);
300                   pmdddldata->SetColumn(icol);
301                   pmdddldata->SetSignal(isig);
302                   pmdddldata->SetParityBit(ibit);
303                   
304                   pmdddlcont->Add(pmdddldata);
305                   
306                 } // data word loop
307
308               if (iwordcount == blRawDataLength) break;
309
310             } // patch bus loop
311
312           if (dspHeader.GetPaddingWord() == 1) iuppLimit++;
313           if (iwordcount == blRawDataLength) break;
314
315         } // end of DSP
316       if (iwordcount == blRawDataLength) break;
317
318     } // end of BLOCK
319
320   
321   delete [] buffer;
322
323   return kTRUE;
324 }
325 //_____________________________________________________________________________
326 void AliPMDRawStream::GetRowCol(Int_t ddlno, Int_t pbusid,
327                                 UInt_t mcmno, UInt_t chno,
328                                 Int_t startRowBus[], Int_t endRowBus[],
329                                 Int_t startColBus[], Int_t endColBus[],
330                                 Int_t &row, Int_t &col) const
331 {
332 // decode: ddlno, patchbusid, mcmno, chno -> um, row, col
333
334
335   static const UInt_t kCh[64] = { 53, 58, 57, 54, 61, 62, 60, 63,
336                                   49, 59, 56, 55, 52, 50, 48, 51,
337                                   44, 47, 45, 43, 40, 39, 36, 46,
338                                   32, 35, 33, 34, 41, 38, 37, 42,
339                                   21, 26, 25, 22, 29, 30, 28, 31,
340                                   17, 27, 24, 23, 20, 18, 16, 19,
341                                   12, 15, 13, 11,  8,  7,  4, 14,
342                                   0,   3,  1,  2,  9,  6,  5, 10 };
343
344
345   Int_t rowcol  = kCh[chno];
346   Int_t irownew = rowcol/4;
347   Int_t icolnew = rowcol%4;
348
349   if (ddlno == 0 )
350     {
351       row = startRowBus[pbusid] + irownew;
352       col = startColBus[pbusid] + mcmno*4 + icolnew;
353     }
354   else if (ddlno == 1)
355     {
356     row = endRowBus[pbusid] - (15 - irownew);
357     col = startColBus[pbusid] + mcmno*4 + icolnew;
358     
359     }
360   else if (ddlno == 2 )
361     {
362       row = startRowBus[pbusid] + irownew;
363       col = endColBus[pbusid] - mcmno*4 - (3 - icolnew);
364     }
365   else if (ddlno == 3)
366     {
367     row = endRowBus[pbusid] - (15 - irownew);
368     col = endColBus[pbusid] - mcmno*4 - (3 - icolnew);
369     }
370   else if (ddlno == 4 )
371     {
372       if (pbusid  < 18)
373         {
374           if(mcmno > 11)
375             {
376               row = startRowBus[pbusid] + irownew;
377               col = startColBus[pbusid] + (mcmno-12)*4 + icolnew;
378             }
379           else
380             {
381               // Add 16 to skip the 1st 15 rows
382               row = startRowBus[pbusid] + irownew + 16;
383               col = startColBus[pbusid] + (mcmno)*4 + icolnew;
384             }
385         }
386       else if(pbusid > 17)
387         {
388           if(mcmno > 11)
389             {
390               row = endRowBus[pbusid] - (15 - irownew)  ;
391               col = endColBus[pbusid] - (mcmno - 12)*4 - (3 - icolnew);
392             }
393           else 
394             {
395               if(endRowBus[pbusid] - startRowBus[pbusid] > 16)
396                 row = endRowBus[pbusid] - (15 - irownew) - 16 ;
397               else
398                 row = endRowBus[pbusid] - (15 - irownew) ;
399               col = endColBus[pbusid] - mcmno*4 - (3 - icolnew); 
400             }
401         }
402     }
403   
404   else if (ddlno == 5)
405     {
406       if (pbusid  <= 17)
407         {
408           if (mcmno > 11)
409             {
410               // Subtract 16 to skip the 1st 15 rows
411               row = endRowBus[pbusid] - 16 -(15 - irownew);
412               col = startColBus[pbusid] + (mcmno-12)*4 + icolnew;
413             }
414           else
415             {
416               row = endRowBus[pbusid]  - (15 - irownew) ;
417               col = startColBus[pbusid] + mcmno*4 + icolnew;
418             }
419         }
420       
421       else if (pbusid > 17)
422         {
423           if(mcmno > 11)
424             {
425               // Add 16 to skip the 1st 15 rows
426               row = startRowBus[pbusid] + irownew + 16;
427               col = endColBus[pbusid] - (mcmno - 12)*4 - (3 - icolnew);
428             }
429           else 
430             {
431               row = startRowBus[pbusid] + irownew ;
432               col = endColBus[pbusid] - mcmno*4 - (3 - icolnew); 
433             }
434         }
435     }
436   
437 }
438 //_____________________________________________________________________________
439 void AliPMDRawStream::ConvertDDL2SMN(Int_t iddl, Int_t imodule,
440                                      Int_t &smn, Int_t &detector) const
441 {
442   // This converts the DDL number (0 to 5), Module Number (0-47)
443   // to Serial module number in one detector (SMN : 0-23) and
444   // detector number (0:PRE plane, 1:CPV plane)
445   if (iddl < 4)
446     {
447       smn = imodule;
448       detector = 0;
449     }
450   else
451     {
452       smn = imodule - 24;
453       detector = 1;
454     }
455 }
456 //_____________________________________________________________________________
457 void AliPMDRawStream::TransformH2S(Int_t smn, Int_t &row, Int_t &col) const
458 {
459   // This does the transformation of the hardware coordinate to
460   // software 
461   // i.e., For SuperModule 0 &1, instead of 96x48(hardware),
462   // it is 48x96 (software)
463   // For Supermodule 3 & 4, 48x96
464
465   Int_t irownew  = 0;
466   Int_t icolnew  = 0;
467
468   if(smn < 12)
469     {
470       irownew = col;
471       icolnew = row;
472     }
473   else if(smn >= 12 && smn < 24)
474     {
475       irownew = row;
476       icolnew = col;
477     }
478
479   row = irownew;
480   col = icolnew;
481 }
482 //_____________________________________________________________________________
483 Int_t AliPMDRawStream::ComputeParity(Int_t data)
484 {
485 // Calculate the parity bit
486
487   Int_t count = 0;
488   for(Int_t j = 0; j<29; j++)
489     {
490       if (data & 0x01 ) count++;
491       data >>= 1;
492     }
493   
494   Int_t parity = count%2;
495
496   return parity;
497 }
498
499 //_____________________________________________________________________________