]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PMD/AliPMDDDLRawData.cxx
Correct number of tracks that is returned by LoadPoints method. Removal of double...
[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 #include <TString.h>
22 #include <TSystem.h>
23
24 #include "AliLog.h"
25 #include "AliRawDataHeaderSim.h"
26 #include "AliBitPacking.h"
27 #include "AliPMDdigit.h"
28 #include "AliPMDBlockHeader.h"
29 #include "AliPMDDspHeader.h"
30 #include "AliPMDPatchBusHeader.h"
31 #include "AliPMDRawStream.h"
32 #include "AliPMDDDLRawData.h"
33 #include "AliDAQ.h"
34 #include "AliFstream.h"
35
36 ClassImp(AliPMDDDLRawData)
37
38 AliPMDDDLRawData::AliPMDDDLRawData():
39   fDigits(new TClonesArray("AliPMDdigit", 1000))
40 {
41   // Default Constructor
42   //
43
44 }
45 //____________________________________________________________________________
46 AliPMDDDLRawData::AliPMDDDLRawData(const AliPMDDDLRawData& ddlraw):
47   TObject(ddlraw),
48   fDigits(ddlraw.fDigits)
49 {
50   //Copy Constructor 
51 }
52 //____________________________________________________________________________
53 AliPMDDDLRawData & AliPMDDDLRawData::operator=(const AliPMDDDLRawData& ddlraw)
54 {
55   //Assignment operator 
56   if(this != &ddlraw)
57     {
58       fDigits = ddlraw.fDigits;
59     }
60   return *this;
61 }
62 //____________________________________________________________________________
63
64 AliPMDDDLRawData::~AliPMDDDLRawData()
65 {
66   // Default Destructor
67   //
68
69 }
70
71 //____________________________________________________________________________
72 void AliPMDDDLRawData::WritePMDRawData(TTree *treeD)
73 {
74   // write digits into raw data format
75
76   AliFstream *outfile;
77
78   TBranch *branch = treeD->GetBranch("PMDDigit");
79   if (!branch)
80     {
81       AliError("PMD Digit branch not found");
82       return;
83     }
84   branch->SetAddress(&fDigits);  
85   
86   Int_t   nmodules = (Int_t) treeD->GetEntries();
87   AliDebug(1,Form("Number of modules inside treeD = %d",nmodules));
88
89   const Int_t kDDL          = AliDAQ::NumberOfDdls("PMD");
90   Int_t modulePerDDL        = 0;
91
92
93   AliRawDataHeaderSim header;
94   UInt_t sizeRawData = 0;
95   
96   const Int_t kbusSize = 51;
97   const Int_t kSize = 1536;
98   UInt_t buffer[kSize];
99
100   UInt_t busPatch[kbusSize][1536];
101
102   Int_t contentsBus[kbusSize];
103
104   Char_t filename[80];
105
106
107   Int_t mmodule = 0;
108   for(Int_t iddl = 0; iddl < kDDL; iddl++)
109     {
110       strcpy(filename,AliDAQ::DdlFileName("PMD",iddl));
111       
112       outfile = new AliFstream(filename);
113       
114       if (iddl < 4)
115         {
116           modulePerDDL = 6;
117           mmodule = 6*iddl;
118         }
119       else if (iddl == 4)
120         {
121           modulePerDDL = 12;
122           mmodule = 24;
123         }
124       else if (iddl == 5)
125         {
126           modulePerDDL = 12;
127           mmodule = 30;
128         }
129
130
131
132       // Write the Dummy Data Header into the file
133       Int_t bHPosition = outfile->Tellp();
134       outfile->WriteBuffer((char*)(&header),sizeof(header));
135
136       for (Int_t ibus = 0; ibus < kbusSize; ibus++)
137         {
138           contentsBus[ibus] = 0;
139           for (Int_t ich = 0; ich < kSize; ich++)
140             {
141               busPatch[ibus][ich] = 0;
142             }
143         }
144
145       for(Int_t ium = 0; ium < modulePerDDL; ium++)
146         {
147           if (iddl == 4 && ium == 6) mmodule = 42;
148
149           // Extract energy deposition per cell and pack it
150           // in a 32-bit word and returns all the total words
151           // per one unit-module
152           
153           GetUMDigitsData(treeD, mmodule, iddl, contentsBus, busPatch);
154           mmodule++;
155         }
156
157       
158
159       Int_t ij = 0;
160       Int_t dsp[10];
161       Int_t dspBus[10];
162       for (Int_t i = 0; i < 10; i++)
163         {
164           dsp[i] = 0;
165           dspBus[i] = 0;
166           for (Int_t ibus=0; ibus < 5; ibus++)
167             {
168               ij++;
169               if (contentsBus[ij] > 0)
170                 {
171                   dsp[i] += contentsBus[ij];
172                   dspBus[i]++;
173                 }
174               //ij++;
175             }
176           // Add the patch Bus header to the DSP contents
177           dsp[i] += 4*dspBus[i];
178         }
179
180       Int_t dspBlockARDL    = 0;
181       Int_t dspBlockBRDL    = 0;
182
183       for (Int_t i = 0; i < 5; i++)
184         {
185           Int_t ieven = 2*i;
186           Int_t iodd  = 2*i + 1;
187           if (dsp[ieven] > 0)
188             {
189               dspBlockARDL += dsp[ieven];
190               Int_t remainder = dsp[ieven]%2;
191               if (remainder == 1)
192                 {
193                   dspBlockARDL++;
194                 }
195             }
196           if (dsp[iodd] > 0)
197             {
198               dspBlockBRDL += dsp[iodd];
199               Int_t remainder = dsp[ieven]%2;
200               if (remainder == 1)
201                 {
202                   dspBlockARDL++;
203                 }
204             }
205         }
206       
207       // Start writing the DDL file
208
209       AliPMDBlockHeader blHeader;
210       AliPMDDspHeader   dspHeader;
211       AliPMDPatchBusHeader pbusHeader;
212
213       const Int_t kblHLen   = blHeader.GetHeaderLength();
214       const Int_t kdspHLen  = dspHeader.GetHeaderLength();
215       const Int_t kpbusHLen = pbusHeader.GetHeaderLength();
216
217       UInt_t dspRDL = 0;
218       UInt_t dspBlockHeaderWord[8];
219       UInt_t dspHeaderWord[10];
220       UInt_t patchBusHeaderWord[4];
221       Int_t  iskip[5];
222
223
224       for (Int_t iblock = 0; iblock < 2; iblock++)
225         {
226           // DSP Block Header
227           
228           for (Int_t i=0; i<kblHLen; i++)
229             {
230               dspBlockHeaderWord[i] = 0;
231             }
232           if (iblock == 0)
233             {
234               dspBlockHeaderWord[1] = (UInt_t) (dspBlockARDL + kblHLen);
235               dspBlockHeaderWord[2] = (UInt_t) dspBlockARDL;
236             }
237           else if (iblock == 1)
238             {
239               dspBlockHeaderWord[1] = (UInt_t) (dspBlockBRDL + kblHLen);
240               dspBlockHeaderWord[2] = (UInt_t) dspBlockBRDL;
241             }
242
243           outfile->WriteBuffer((char*)dspBlockHeaderWord,kblHLen*sizeof(UInt_t));
244
245           if (iblock == 0)
246             {
247               iskip[0] = 0;
248               iskip[1] = 10;
249               iskip[2] = 20;
250               iskip[3] = 30;
251               iskip[4] = 40;
252             }
253           else if (iblock == 1)
254             {
255               iskip[0] = 5;
256               iskip[1] = 15;
257               iskip[2] = 25;
258               iskip[3] = 35;
259               iskip[4] = 45;
260             }
261
262           for (Int_t idsp = 0; idsp < 5; idsp++)
263             {
264               // DSP Header
265               Int_t dspno = 0;
266               if (iblock == 0)
267                 {
268                   dspno = 2*idsp;
269                   dspRDL = (UInt_t) dsp[dspno];
270                 }
271               else if (iblock == 1)
272                 {
273                   dspno = 2*idsp + 1;
274                   dspRDL = (UInt_t) dsp[dspno];
275                 }
276
277               for (Int_t i=0; i<kdspHLen; i++)
278                 {
279                   dspHeaderWord[i] = 0;
280                 }
281               Int_t remainder = dspRDL%2;
282               if (remainder == 1) dspRDL++;
283
284               dspHeaderWord[1] = dspRDL + kdspHLen;
285               dspHeaderWord[2] = dspRDL;
286               dspHeaderWord[3] = dspno;
287               if (remainder == 1) dspHeaderWord[8] = 1; // setting the padding word
288               outfile->WriteBuffer((char*)dspHeaderWord,kdspHLen*sizeof(UInt_t));
289
290               for (Int_t ibus = 0; ibus < 5; ibus++)
291                 {
292                   // Patch Bus Header
293                   // BKN - just added 1
294                   Int_t busno = iskip[idsp] + ibus + 1;
295                   Int_t patchbusRDL = contentsBus[busno];
296
297                   if (patchbusRDL > 0)
298                     {
299                       patchBusHeaderWord[0] = 0;
300                       patchBusHeaderWord[1] = (UInt_t) (patchbusRDL + kpbusHLen);
301                       patchBusHeaderWord[2] = (UInt_t) patchbusRDL;
302                       patchBusHeaderWord[3] = (UInt_t) busno;
303                     }
304                   else if (patchbusRDL == 0)
305                     {
306                       patchBusHeaderWord[0] = 0;
307                       patchBusHeaderWord[1] = (UInt_t) kpbusHLen;
308                       patchBusHeaderWord[2] = (UInt_t) 0;
309                       patchBusHeaderWord[3] = (UInt_t) busno;
310                     }
311
312
313                   outfile->WriteBuffer((char*)patchBusHeaderWord,4*sizeof(UInt_t));
314
315
316                   for (Int_t iword = 0; iword < patchbusRDL; iword++)
317                     {
318                       buffer[iword] = busPatch[busno][iword];
319                     }
320                   
321                   outfile->WriteBuffer((char*)buffer,patchbusRDL*sizeof(UInt_t));
322
323                 } // End of patch bus loop
324
325
326               // Adding a padding word if the total words odd
327               if (remainder == 1)
328                 {
329                   UInt_t paddingWord = dspHeader.GetDefaultPaddingWord();
330                   outfile->WriteBuffer((char*)(&paddingWord),sizeof(UInt_t));
331                 }
332             }
333         }
334
335       // Write real data header
336       // take the pointer to the beginning of the data header
337       // write the total number of words per ddl and bring the
338       // pointer to the current file position and close it
339       UInt_t cFPosition = outfile->Tellp();
340       sizeRawData = cFPosition - bHPosition - sizeof(header);
341
342       header.fSize = cFPosition - bHPosition;
343       header.SetAttribute(0);  // valid data
344       outfile->Seekp(bHPosition);
345       outfile->WriteBuffer((char*)(&header),sizeof(header));
346       outfile->Seekp(cFPosition);
347
348       delete outfile;
349     } // DDL Loop over
350
351
352 }
353 //____________________________________________________________________________
354 void AliPMDDDLRawData::GetUMDigitsData(TTree *treeD, Int_t imodule,
355                                        Int_t ddlno, Int_t *contentsBus,
356                                        UInt_t busPatch[][1536])
357 {
358   // Retrives digits data UnitModule by UnitModule
359
360   UInt_t baseword;
361   UInt_t mcmno, chno;
362   UInt_t adc;
363   Int_t  det, smn, irow, icol;
364   Int_t  parity;
365   
366   const Int_t kMaxBus = 51;   // BKN
367   Int_t totPatchBus, bPatchBus, ePatchBus;
368   Int_t ibus, totmcm, rows, cols, rowe, cole;
369   Int_t moduleno;
370   Int_t busno = 0;
371   Int_t patchBusNo[kMaxBus], mcmperBus[kMaxBus];
372   Int_t startRowBus[kMaxBus], startColBus[kMaxBus];
373   Int_t endRowBus[kMaxBus], endColBus[kMaxBus];
374
375   Int_t beginPatchBus = -1;
376   Int_t endPatchBus   = -1;
377   for(Int_t i = 0; i < kMaxBus; i++)
378     {
379       patchBusNo[i]  = -1;
380       mcmperBus[i]   = -1;
381       startRowBus[i] = -1;
382       startColBus[i] = -1;
383       endRowBus[i]   = -1;
384       endColBus[i]   = -1;
385     }
386   Int_t modulePerDDL = 0;
387   if (ddlno < 4)
388     {
389       modulePerDDL = 6;
390     }
391   else if (ddlno == 4 || ddlno == 5)
392     {
393       modulePerDDL = 12;
394     }
395
396   TString fileName(gSystem->Getenv("ALICE_ROOT"));
397
398   if(ddlno == 0)
399     {
400       fileName += "/PMD/PMD_Mapping_ddl0.dat";
401     }
402   else if(ddlno == 1)
403     {
404       fileName += "/PMD/PMD_Mapping_ddl1.dat";
405     }
406   else if(ddlno == 2)
407     {
408       fileName += "/PMD/PMD_Mapping_ddl2.dat";
409     }
410   else if(ddlno == 3)
411     {
412       fileName += "/PMD/PMD_Mapping_ddl3.dat";
413     }
414   else if(ddlno == 4)
415     {
416       fileName += "/PMD/PMD_Mapping_ddl4.dat";
417     }
418   else if(ddlno == 5)
419     {
420       fileName += "/PMD/PMD_Mapping_ddl5.dat";
421     }
422
423   ifstream infile;
424   infile.open(fileName.Data(), ios::in); // ascii file
425   if(!infile)
426     AliError(Form("Could not read the mapping file for DDL No = %d",ddlno));
427
428   for (Int_t im = 0; im < modulePerDDL; im++)
429     {
430       infile >> moduleno;
431       infile >> totPatchBus >> bPatchBus >> ePatchBus;
432
433       if (moduleno == imodule)
434         {
435           beginPatchBus = bPatchBus;
436           endPatchBus   = ePatchBus;
437         }
438       
439       for(Int_t i=0; i<totPatchBus; i++)
440         {
441           infile >> ibus >> totmcm >> rows >> rowe >> cols >> cole;
442
443           if (moduleno == imodule)
444             {
445               patchBusNo[ibus]   = ibus;
446               mcmperBus[ibus]    = totmcm;
447               startRowBus[ibus]  = rows;
448               startColBus[ibus]  = cols;
449               endRowBus[ibus]    = rowe;
450               endColBus[ibus]    = cole;
451             }
452         }
453
454     }
455
456   infile.close();
457
458   treeD->GetEntry(imodule); 
459   Int_t nentries = fDigits->GetLast();
460   Int_t totword = nentries+1;
461
462   AliPMDdigit *pmddigit;
463
464   for (Int_t ient = 0; ient < totword; ient++)
465     {
466       pmddigit = (AliPMDdigit*)fDigits->UncheckedAt(ient);
467       
468       det    = pmddigit->GetDetector();
469       smn    = pmddigit->GetSMNumber();
470       irow   = pmddigit->GetRow();
471       icol   = pmddigit->GetColumn();
472       adc    = (UInt_t) pmddigit->GetADC();
473
474       TransformS2H(smn,irow,icol);
475
476       GetMCMCh(ddlno, smn, irow, icol, beginPatchBus, endPatchBus,
477                mcmperBus, startRowBus, startColBus,
478                endRowBus, endColBus, busno, mcmno, chno);
479
480       baseword = 0;
481       AliBitPacking::PackWord(adc,baseword,0,11);
482       AliBitPacking::PackWord(chno,baseword,12,17);
483       AliBitPacking::PackWord(mcmno,baseword,18,28);
484       AliBitPacking::PackWord(0,baseword,29,30);
485       parity = ComputeParity(baseword);      // generate the parity bit
486       AliBitPacking::PackWord(parity,baseword,31,31);
487
488       Int_t jj = contentsBus[busno];
489       busPatch[busno][jj] = baseword;
490
491       contentsBus[busno]++;
492     }
493
494 }
495
496 //____________________________________________________________________________
497 void AliPMDDDLRawData::TransformS2H(Int_t smn, Int_t &irow, Int_t &icol)
498 {
499   // Does the Software to Hardware coordinate transformation
500   //
501
502   Int_t  irownew = 0;
503   Int_t  icolnew = 0;
504
505   // First in digits we have all dimension 48x96
506   // Transform into the realistic one, i.e, For SM 0&1 96(row)x48(col)
507   // and for SM 2&3 48(row)x96(col)
508   // 
509   if(smn < 12)
510     {
511       irownew = icol;
512       icolnew = irow;
513     }
514   else if( smn >= 12 && smn < 24)
515     {
516       irownew = irow;
517       icolnew = icol;
518     }
519
520   irow = irownew;
521   icol = icolnew;
522
523 }
524
525
526 //____________________________________________________________________________
527
528 void AliPMDDDLRawData::GetMCMCh(Int_t ddlno, Int_t smn, Int_t row, Int_t col,
529                                 Int_t beginPatchBus, Int_t endPatchBus,
530                                 Int_t *mcmperBus,
531                                 Int_t *startRowBus, Int_t *startColBus,
532                                 Int_t *endRowBus, Int_t *endColBus,
533                                 Int_t & busno, UInt_t &mcmno, UInt_t &chno)
534 {
535   // This converts row col to hardware channel number
536   // This is the final version of mapping supplied by Mriganka
537
538     UInt_t iCh[16][4];
539
540 //    ChMap(ddlno, smn, iCh);
541
542
543     static const UInt_t kChDdl01[16][4] = { {6, 4, 5, 7},
544                                            {10, 2, 1, 9},
545                                            {12, 0, 3, 11},
546                                            {14, 8, 13, 15},
547                                            {16, 18, 23, 17},
548                                            {20, 28, 31, 19},
549                                            {22, 30, 29, 21},
550                                            {24, 26, 27, 25},
551                                            {38, 36, 37, 39},
552                                            {42, 34, 33, 41},
553                                            {44, 32, 35, 43},
554                                            {46, 40, 45, 47},
555                                            {48, 50, 55, 49},
556                                            {52, 60, 63, 51},
557                                            {54, 62, 61, 53},
558                                            {56, 58, 59, 57} };
559
560
561     static const UInt_t kChDdl23[16][4] = { {57, 59, 58, 56},
562                                             {53, 61, 62, 54},
563                                             {51, 63, 60, 52},
564                                             {49, 55, 50, 48},
565                                             {47, 45, 40, 46},
566                                             {43, 35, 32, 44},
567                                             {41, 33, 34, 42},
568                                             {39, 37, 36, 38},
569                                             {25, 27, 26, 24},
570                                             {21, 29, 30, 22},
571                                             {19, 31, 28, 20},
572                                             {17, 23, 18, 16},
573                                             {15, 13, 8, 14},
574                                             {11, 3, 0, 12},
575                                             {9, 1, 2, 10},
576                                             {7, 5, 4, 6} };
577     
578     
579     static const UInt_t kChDdl41[16][4] = { {56, 58, 59, 57},
580                                            {54, 62, 61, 53},
581                                            {52, 60, 63, 51},
582                                            {48, 50, 55, 49},
583                                            {46, 40, 45, 47},
584                                            {44, 32, 35, 43},
585                                            {42, 34, 33, 41},
586                                            {38, 36, 37, 39},
587                                            {24, 26, 27, 25},
588                                            {22, 30, 29, 21},
589                                            {20, 28, 31, 19},
590                                            {16, 18, 23, 17},
591                                            {14, 8, 13, 15},
592                                            {12, 0, 3, 11},
593                                            {10, 2, 1, 9},
594                                            {6, 4, 5, 7} };
595
596
597     static const UInt_t kChDdl42[16][4] = { {7, 5, 4, 6},
598                                             {9, 1, 2, 10},
599                                             {11, 3, 0, 12},
600                                             {15, 13, 8, 14},
601                                             {17, 23, 18, 16},
602                                             {19, 31, 28, 20},
603                                             {21, 29, 30, 22},
604                                             {25, 27, 26, 24},
605                                             {39, 37, 36, 38},
606                                             {41, 33, 34, 42},
607                                             {43, 35, 32, 44},
608                                             {47, 45, 40, 46},
609                                             {49, 55, 50, 48},
610                                             {51, 63, 60, 52},
611                                             {53, 61, 62, 54},
612                                             {57, 59, 58, 56} };
613
614
615     static const UInt_t kChDdl51[16][4] = { {7, 5, 4, 6},
616                                             {9, 1, 2, 10},
617                                             {11, 3, 0, 12},
618                                             {15, 13, 8, 14},
619                                             {17, 23, 18, 16},
620                                             {19, 31, 28, 20},
621                                             {21, 29, 30, 22},
622                                             {25, 27, 26, 24},
623                                             {39, 37, 36, 38},
624                                             {41, 33, 34, 42},
625                                             {43, 35, 32, 44},
626                                             {47, 45, 40, 46},
627                                             {49, 55, 50, 48},
628                                             {51, 63, 60, 52},
629                                             {53, 61, 62, 54},
630                                             {57, 59, 58, 56} };
631     
632
633
634     static const UInt_t kChDdl52[16][4] = { {56, 58, 59, 57},
635                                             {54, 62, 61, 53},
636                                             {52, 60, 63, 51},
637                                             {48, 50, 55, 49},
638                                             {46, 40, 45, 47},
639                                             {44, 32, 35, 43},
640                                             {42, 34, 33, 41},
641                                             {38, 36, 37, 39},
642                                             {24, 26, 27, 25},
643                                             {22, 30, 29, 21},
644                                             {20, 28, 31, 19},
645                                             {16, 18, 23, 17},
646                                             {14, 8, 13, 15},
647                                             {12, 0, 3, 11},
648                                             {10, 2, 1, 9},
649                                             {6, 4, 5, 7} };
650     
651     
652     for (Int_t i = 0; i < 16; i++)
653     {
654         for (Int_t j = 0; j < 4; j++)
655         {
656             if (ddlno == 0 || ddlno == 1) iCh[i][j] = kChDdl01[i][j];
657             if (ddlno == 2 || ddlno == 3) iCh[i][j] = kChDdl23[i][j];
658             
659             if (ddlno == 4 && smn < 6)                iCh[i][j] = kChDdl41[i][j];
660             if (ddlno == 4 && (smn >= 18 && smn < 24))iCh[i][j] = kChDdl42[i][j];
661             if (ddlno == 5 && (smn >= 12 && smn < 18))iCh[i][j] = kChDdl51[i][j];
662             if (ddlno == 5 && (smn >=  6 && smn < 12))iCh[i][j] = kChDdl52[i][j];
663         }
664     }
665
666
667   Int_t irownew = row%16;
668   Int_t icolnew = col%4;
669   
670   chno  = iCh[irownew][icolnew];
671   
672   
673   for (Int_t ibus = beginPatchBus; ibus <= endPatchBus; ibus++)
674     {
675       Int_t srow = startRowBus[ibus];
676       Int_t erow = endRowBus[ibus];
677       Int_t scol = startColBus[ibus];
678       Int_t ecol = endColBus[ibus];
679       Int_t tmcm = mcmperBus[ibus];
680       if ((row >= srow && row <= erow) && (col >= scol && col <= ecol))
681         {
682           busno = ibus;
683           
684           // Find out the MCM Number
685           //
686           
687           if (ddlno == 0 || ddlno == 1)
688             {
689               // PRE plane, SU Mod = 0, 1
690               mcmno = (col-scol)/4 + 1;
691               
692             }
693           else if (ddlno == 2 || ddlno == 3)
694             {
695               // PRE plane,  SU Mod = 2, 3
696               Int_t icolnew = (col - scol)/4;
697               mcmno = tmcm - icolnew;
698             }
699           else if (ddlno == 4 )
700             {
701               // CPV plane,  SU Mod = 0, 3 : ddl = 4
702               
703               if(ibus <= 18)
704                 {
705                   Int_t midrow = srow + 16;
706                   if(row >= srow && row < midrow)
707                     {
708                       mcmno = 12 + (col-scol)/4 + 1;
709                     }
710                   else if(row >= midrow && row <= erow)
711                   
712                     {
713                       mcmno = (col-scol)/4 + 1;
714                     }
715                 }
716               
717               else if (ibus > 18)
718                 {
719                   Int_t rowdiff = endRowBus[ibus] - startRowBus[ibus];
720                   if(rowdiff > 16)
721                     {
722                       Int_t midrow = srow + 16;
723                       if (row >= midrow && row <= erow)
724                         {
725                           mcmno = 12 + (ecol -col)/4 + 1;
726                         }
727                       else if (row >= srow && row < midrow)
728                         {
729                           mcmno = (ecol - col)/4 + 1;
730                         }
731                     }
732                   else if (rowdiff < 16)
733                     {
734                       mcmno = (ecol - col)/4 + 1;
735                     }
736                 }
737             }
738           else if ( ddlno == 5)
739             {
740               // CPV plane,  SU Mod = 1, 2 : ddl = 5
741               
742               if(ibus <= 18)
743                 {
744                   Int_t midrow = srow + 16;
745                   if(row >= srow && row < midrow)
746                     {
747                       mcmno = 12 + (col-scol)/4 + 1;
748                     }
749                   else if(row >= midrow && row <= erow)
750                     {
751                       mcmno = (col-scol)/4 + 1;
752                     }
753                 }
754               
755               else if (ibus > 18)
756                 {
757                   Int_t rowdiff = endRowBus[ibus] - startRowBus[ibus];
758                   if(rowdiff > 16)
759                     {
760                       Int_t midrow = srow + 16;
761                       if (row >= midrow && row <= erow)
762                         {
763                           mcmno = 12 + (ecol - col)/4 + 1;
764                         }
765                       else if (row >= srow && row < midrow)
766                         {
767                           mcmno = (ecol - col)/4 + 1;
768                         }
769                     }
770                   else if (rowdiff < 16)
771                     {
772                       mcmno = (ecol - col)/4 + 1;
773                     }
774                 }
775             }
776         }
777     }
778
779
780 //____________________________________________________________________________
781
782 /*
783 void ChMap(Int_t ddlno, Int_t smn, UInt_t iCh[][4])
784 {
785 // Row and Col converted to Channel
786
787     static const UInt_t kChDdl01[16][4] = { {6, 4, 5, 7},
788                                            {10, 2, 1, 9},
789                                            {12, 0, 3, 11},
790                                            {14, 8, 13, 15},
791                                            {16, 18, 23, 17},
792                                            {20, 28, 31, 19},
793                                            {22, 30, 29, 21},
794                                            {24, 26, 27, 25},
795                                            {38, 36, 37, 39},
796                                            {42, 34, 33, 41},
797                                            {44, 32, 35, 43},
798                                            {46, 40, 45, 47},
799                                            {48, 50, 55, 49},
800                                            {52, 60, 63, 51},
801                                            {54, 62, 61, 53},
802                                            {56, 58, 59, 57} };
803
804
805     static const UInt_t kChDdl23[16][4] = { {57, 59, 58, 56},
806                                             {53, 61, 62, 54},
807                                             {51, 63, 60, 52},
808                                             {49, 55, 50, 48},
809                                             {47, 45, 40, 46},
810                                             {43, 35, 32, 44},
811                                             {41, 33, 34, 42},
812                                             {39, 37, 36, 38},
813                                             {25, 27, 26, 24},
814                                             {21, 29, 30, 22},
815                                             {19, 31, 28, 20},
816                                             {17, 23, 18, 16},
817                                             {15, 13, 8, 14},
818                                             {11, 3, 0, 12},
819                                             {9, 1, 2, 10},
820                                             {7, 5, 4, 6} };
821     
822     
823     static const UInt_t kChDdl41[16][4] = { {56, 58, 59, 57},
824                                            {54, 62, 61, 53},
825                                            {52, 60, 63, 51},
826                                            {48, 50, 55, 49},
827                                            {46, 40, 45, 47},
828                                            {44, 32, 35, 43},
829                                            {42, 34, 33, 41},
830                                            {38, 36, 37, 39},
831                                            {24, 26, 27, 25},
832                                            {22, 30, 29, 21},
833                                            {20, 28, 31, 19},
834                                            {16, 18, 23, 17},
835                                            {14, 8, 13, 15},
836                                            {12, 0, 3, 11},
837                                            {10, 2, 1, 9},
838                                            {6, 4, 5, 7} };
839
840
841     static const UInt_t kChDdl42[16][4] = { {7, 5, 4, 6},
842                                             {9, 1, 2, 10},
843                                             {11, 3, 0, 12},
844                                             {15, 13, 8, 14},
845                                             {17, 23, 18, 16},
846                                             {19, 31, 28, 20},
847                                             {21, 29, 30, 22},
848                                             {25, 27, 26, 24},
849                                             {39, 37, 36, 38},
850                                             {41, 33, 34, 42},
851                                             {43, 35, 32, 44},
852                                             {47, 45, 40, 46},
853                                             {49, 55, 50, 48},
854                                             {51, 63, 60, 52},
855                                             {53, 61, 62, 54},
856                                             {57, 59, 58, 56} };
857
858
859     static const UInt_t kChDdl51[16][4] = { {7, 5, 4, 6},
860                                             {9, 1, 2, 10},
861                                             {11, 3, 0, 12},
862                                             {15, 13, 8, 14},
863                                             {17, 23, 18, 16},
864                                             {19, 31, 28, 20},
865                                             {21, 29, 30, 22},
866                                             {25, 27, 26, 24},
867                                             {39, 37, 36, 38},
868                                             {41, 33, 34, 42},
869                                             {43, 35, 32, 44},
870                                             {47, 45, 40, 46},
871                                             {49, 55, 50, 48},
872                                             {51, 63, 60, 52},
873                                             {53, 61, 62, 54},
874                                             {57, 59, 58, 56} };
875     
876
877
878     static const UInt_t kChDdl52[16][4] = { {56, 58, 59, 57},
879                                             {54, 62, 61, 53},
880                                             {52, 60, 63, 51},
881                                             {48, 50, 55, 49},
882                                             {46, 40, 45, 47},
883                                             {44, 32, 35, 43},
884                                             {42, 34, 33, 41},
885                                             {38, 36, 37, 39},
886                                             {24, 26, 27, 25},
887                                             {22, 30, 29, 21},
888                                             {20, 28, 31, 19},
889                                             {16, 18, 23, 17},
890                                             {14, 8, 13, 15},
891                                             {12, 0, 3, 11},
892                                             {10, 2, 1, 9},
893                                             {6, 4, 5, 7} };
894     
895     
896     for (Int_t i = 0; i < 16; i++)
897     {
898         for (Int_t j = 0; j < 4; j++)
899         {
900             if (ddlno == 0 || ddlno == 1) iCh[i][j] = kChDdl01[i][j];
901             if (ddlno == 2 || ddlno == 3) iCh[i][j] = kChDdl23[i][j];
902             
903             if (ddlno == 4 && smn < 6)                iCh[i][j] = kChDdl41[i][j];
904             if (ddlno == 4 && (smn >= 18 && smn < 24))iCh[i][j] = kChDdl42[i][j];
905             if (ddlno == 5 && (smn >= 12 && smn < 18))iCh[i][j] = kChDdl51[i][j];
906             if (ddlno == 0 && (smn >=  6 && smn < 12))iCh[i][j] = kChDdl52[i][j];
907         }
908     }
909     
910 }
911 */
912 //____________________________________________________________________________
913
914
915
916 Int_t AliPMDDDLRawData::ComputeParity(UInt_t baseword)
917 {
918   // Generate the parity bit
919
920   Int_t count = 0;
921   for(Int_t j=0; j<29; j++)
922     {
923       if (baseword & 0x01 ) count++;
924       baseword >>= 1;
925     }
926   Int_t parity = count%2;
927   return parity;
928 }
929
930 //____________________________________________________________________________