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