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