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