]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PMD/AliPMDRawStream.cxx
bugs fixed
[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   for (Int_t ibus = 0; ibus < kNPatchBus; ibus++)
164     {
165       mcmperBus[ibus]   = -1;
166       startRowBus[ibus] = -1;
167       endRowBus[ibus]   = -1;
168       startColBus[ibus] = -1;
169       endColBus[ibus]   = -1;
170     }
171
172   for (Int_t im = 0; im < modulePerDDL; im++)
173     {
174       infile >> modno;
175       infile >> totPatchBus >> bPatchBus >> ePatchBus;
176       
177       for(Int_t i=0; i<totPatchBus; i++)
178         {
179           infile >> ibus >> totmcm >> rows >> rowe >> cols >> cole;
180
181           moduleNo[ibus]     = modno;
182           mcmperBus[ibus]    = totmcm;
183           startRowBus[ibus]  = rows;
184           endRowBus[ibus]    = rowe;
185           startColBus[ibus]  = cols;
186           endColBus[ibus]    = cole;
187         }
188     }
189
190   infile.close();
191
192   AliPMDBlockHeader    blockHeader;
193   AliPMDDspHeader      dspHeader;
194   AliPMDPatchBusHeader pbusHeader;
195
196   const Int_t kblHLen   = blockHeader.GetHeaderLength();
197   const Int_t kdspHLen  = dspHeader.GetHeaderLength();
198   const Int_t kpbusHLen = pbusHeader.GetHeaderLength();
199   
200   Int_t parity;
201   Int_t idet, ismn;
202   Int_t irow = -1;
203   Int_t icol = -1;
204
205   Int_t blHeaderWord[8];
206   Int_t dspHeaderWord[10];
207   Int_t pbusHeaderWord[4];
208
209   Int_t ilowLimit       = 0;
210   Int_t iuppLimit       = 0;
211   Int_t blRawDataLength = 0;
212   Int_t iwordcount      = 0;
213
214
215   for (Int_t iblock = 0; iblock < 2; iblock++)
216     {
217       ilowLimit = iuppLimit;
218       iuppLimit = ilowLimit + kblHLen;
219
220       for (Int_t i = ilowLimit; i < iuppLimit; i++)
221         {
222           blHeaderWord[i-ilowLimit] = (Int_t) buffer[i];
223         }
224
225       blockHeader.SetHeader(blHeaderWord);
226
227       blRawDataLength = blockHeader.GetRawDataLength();
228
229       for (Int_t idsp = 0; idsp < 5; idsp++)
230         {
231           ilowLimit = iuppLimit;
232           iuppLimit = ilowLimit + kdspHLen;
233
234           for (Int_t i = ilowLimit; i < iuppLimit; i++)
235             {
236               iwordcount++;
237               dspHeaderWord[i-ilowLimit] = (Int_t) buffer[i];
238             }
239           dspHeader.SetHeader(dspHeaderWord);
240
241           for (Int_t ibus = 0; ibus < 5; ibus++)
242             {
243               ilowLimit = iuppLimit;
244               iuppLimit = ilowLimit + kpbusHLen;
245
246               for (Int_t i = ilowLimit; i < iuppLimit; i++)
247                 {
248                   iwordcount++;
249                   pbusHeaderWord[i-ilowLimit] = (Int_t) buffer[i];
250                 }
251               pbusHeader.SetHeader(pbusHeaderWord);
252               Int_t rawdatalength = pbusHeader.GetRawDataLength();
253               Int_t pbusid = pbusHeader.GetPatchBusId();
254
255               ilowLimit = iuppLimit;
256               iuppLimit = ilowLimit + rawdatalength;
257
258               Int_t imodule = moduleNo[pbusid];
259
260
261               for (Int_t iword = ilowLimit; iword < iuppLimit; iword++)
262                 {
263                   iwordcount++;
264                   data = buffer[iword];
265
266                   Int_t isig =  data & 0x0FFF;
267                   Int_t ich  = (data >> 12) & 0x003F;
268                   Int_t imcm = (data >> 18) & 0x07FF;
269                   Int_t ibit = (data >> 31) & 0x0001;
270                   parity = ComputeParity(data);
271                   if (ibit != parity)
272                     {
273                       AliWarning("ComputeParity:: Parity Error");
274                     }
275                   GetRowCol(iddl, pbusid, imcm, ich, 
276                             startRowBus, endRowBus,
277                             startColBus, endColBus,
278                             irow, icol);
279
280                   ConvertDDL2SMN(iddl, imodule, ismn, idet);
281                   TransformH2S(ismn, irow, icol);
282
283                   pmdddldata = new AliPMDddldata();
284
285                   pmdddldata->SetDetector(idet);
286                   pmdddldata->SetSMN(ismn);
287                   pmdddldata->SetModule(imodule);
288                   pmdddldata->SetPatchBusId(pbusid);
289                   pmdddldata->SetMCM(imcm);
290                   pmdddldata->SetChannel(ich);
291                   pmdddldata->SetRow(irow);
292                   pmdddldata->SetColumn(icol);
293                   pmdddldata->SetSignal(isig);
294                   pmdddldata->SetParityBit(ibit);
295                   
296                   pmdddlcont->Add(pmdddldata);
297                   
298                 } // data word loop
299
300               if (iwordcount == blRawDataLength) break;
301
302             } // patch bus loop
303
304           if (dspHeader.GetPaddingWord() == 1) iuppLimit++;
305           if (iwordcount == blRawDataLength) break;
306
307         } // end of DSP
308       if (iwordcount == blRawDataLength) break;
309
310     } // end of BLOCK
311   
312   delete [] buffer;
313
314   return kTRUE;
315 }
316 //_____________________________________________________________________________
317 void AliPMDRawStream::GetRowCol(Int_t ddlno, Int_t pbusid,
318                                 UInt_t mcmno, UInt_t chno,
319                                 Int_t startRowBus[], Int_t endRowBus[],
320                                 Int_t startColBus[], Int_t endColBus[],
321                                 Int_t &row, Int_t &col) const
322 {
323 // decode: ddlno, patchbusid, mcmno, chno -> um, row, col
324
325
326   static const UInt_t kCh[64] = { 53, 58, 57, 54, 61, 62, 60, 63,
327                                   49, 59, 56, 55, 52, 50, 48, 51,
328                                   44, 47, 45, 43, 40, 39, 36, 46,
329                                   32, 35, 33, 34, 41, 38, 37, 42,
330                                   21, 26, 25, 22, 29, 30, 28, 31,
331                                   17, 27, 24, 23, 20, 18, 16, 19,
332                                   12, 15, 13, 11,  8,  7,  4, 14,
333                                   0,   3,  1,  2,  9,  6,  5, 10 };
334
335
336   Int_t rowcol  = kCh[chno];
337   Int_t irownew = rowcol/4;
338   Int_t icolnew = rowcol%4;
339
340   if (ddlno == 0 )
341     {
342       row = startRowBus[pbusid] + irownew;
343       col = startColBus[pbusid] + mcmno*4 + icolnew;
344     }
345   else if (ddlno == 1)
346     {
347     row = endRowBus[pbusid] - (15 - irownew);
348     col = startColBus[pbusid] + mcmno*4 + icolnew;
349     
350     }
351   else if (ddlno == 2 )
352     {
353       row = startRowBus[pbusid] + irownew;
354       col = endColBus[pbusid] - mcmno*4 - (3 - icolnew);
355     }
356   else if (ddlno == 3)
357     {
358     row = endRowBus[pbusid] - (15 - irownew);
359     col = endColBus[pbusid] - mcmno*4 - (3 - icolnew);
360     }
361   else if (ddlno == 4 )
362     {
363       if (pbusid  < 18)
364         {
365           if (mcmno <= 11)
366             {
367               // Add 16 to skip the 1st 15 rows
368               row = startRowBus[pbusid] + irownew + 16;
369               col = startColBus[pbusid] + (mcmno)*4 + icolnew;
370             }
371           else if(mcmno > 11)
372             {
373               row = startRowBus[pbusid] + irownew;
374               col = startColBus[pbusid] + (mcmno-12)*4 + icolnew;
375             }
376         }
377       else if(pbusid > 17)
378         {
379           if (mcmno <= 11)
380             {
381               col = endColBus[pbusid] - mcmno*4 - (3 - icolnew); 
382
383               if(endRowBus[pbusid] - startRowBus[pbusid] > 16)
384                 row = endRowBus[pbusid] - (15 - irownew) - 16 ;
385               else
386                 row = endRowBus[pbusid] - (15 - irownew) ;
387               
388
389             }
390           else if(mcmno > 11)
391             {
392               row = endRowBus[pbusid] - (15 - irownew)  ;
393               col = endColBus[pbusid] - (mcmno - 12)*4 - (3 - icolnew);
394             }
395         }
396     }
397   else if (ddlno == 5)
398     {
399       if (pbusid  <= 17)
400         {
401           if (mcmno > 11)
402             {
403               // Subtract 16 to skip the 1st 15 rows
404               row = endRowBus[pbusid] - 16 -(15 - irownew);
405               col = startColBus[pbusid] + (mcmno-12)*4 + icolnew;
406             }
407           else
408             {
409               row = endRowBus[pbusid]  - (15 - irownew) ;
410               col = startColBus[pbusid] + mcmno*4 + icolnew;
411             }
412         }
413       
414       else if (pbusid > 17)
415         {
416           if(mcmno > 11)
417             {
418               // Add 16 to skip the 1st 15 rows
419               row = startRowBus[pbusid] + irownew + 16;
420               col = endColBus[pbusid] - (mcmno - 12)*4 - (3 - icolnew);
421             }
422           else 
423             {
424               row = startRowBus[pbusid] + irownew ;
425               col = endColBus[pbusid] - mcmno*4 - (3 - icolnew); 
426             }
427         }
428     }
429   
430 }
431 //_____________________________________________________________________________
432 void AliPMDRawStream::ConvertDDL2SMN(Int_t iddl, Int_t imodule,
433                                      Int_t &smn, Int_t &detector) const
434 {
435   // This converts the DDL number (0 to 5), Module Number (0-47)
436   // to Serial module number in one detector (SMN : 0-23) and
437   // detector number (0:PRE plane, 1:CPV plane)
438   if (iddl < 4)
439     {
440       smn = imodule;
441       detector = 0;
442     }
443   else
444     {
445       smn = imodule - 24;
446       detector = 1;
447     }
448 }
449 //_____________________________________________________________________________
450 void AliPMDRawStream::TransformH2S(Int_t smn, Int_t &row, Int_t &col) const
451 {
452   // This does the transformation of the hardware coordinate to
453   // software 
454   // i.e., For SuperModule 0 &1, instead of 96x48(hardware),
455   // it is 48x96 (software)
456   // For Supermodule 3 & 4, 48x96
457
458   Int_t irownew  = 0;
459   Int_t icolnew  = 0;
460
461   if(smn < 12)
462     {
463       irownew = col;
464       icolnew = row;
465     }
466   else if(smn >= 12 && smn < 24)
467     {
468       irownew = row;
469       icolnew = col;
470     }
471
472   row = irownew;
473   col = icolnew;
474 }
475 //_____________________________________________________________________________
476 Int_t AliPMDRawStream::ComputeParity(Int_t data)
477 {
478 // Calculate the parity bit
479
480   Int_t count = 0;
481   for(Int_t j = 0; j<29; j++)
482     {
483       if (data & 0x01 ) count++;
484       data >>= 1;
485     }
486   
487   Int_t parity = count%2;
488
489   return parity;
490 }
491
492 //_____________________________________________________________________________