]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PMD/AliPMDDDLRawData.cxx
- Rebbined ratios in the TOF range
[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 "AliPMDddlinfoData.h"
33 #include "AliPMDMappingData.h"
34 #include "AliPMDDDLRawData.h"
35 #include "AliDAQ.h"
36 #include "AliFstream.h"
37 #include "AliCDBManager.h"
38 #include "AliCDBStorage.h"
39 #include "AliCDBEntry.h"
40
41 ClassImp(AliPMDDDLRawData)
42
43 AliPMDDDLRawData::AliPMDDDLRawData():
44   fDdlinfo(GetDdlinfoData()),
45   fMapData(GetMappingData()),
46   fDigits(new TClonesArray("AliPMDdigit", 1000))
47 {
48   // Default Constructor
49   //
50
51 }
52 //____________________________________________________________________________
53 AliPMDDDLRawData::AliPMDDDLRawData(const AliPMDDDLRawData& ddlraw):
54   TObject(ddlraw),
55   fDdlinfo(ddlraw.fDdlinfo),
56   fMapData(ddlraw.fMapData),
57   fDigits(ddlraw.fDigits)
58 {
59   //Copy Constructor 
60 }
61 //____________________________________________________________________________
62 AliPMDDDLRawData & AliPMDDDLRawData::operator=(const AliPMDDDLRawData& ddlraw)
63 {
64   //Assignment operator 
65   if(this != &ddlraw)
66     {
67       fDdlinfo = ddlraw.fDdlinfo;
68       fMapData = ddlraw.fMapData;
69       fDigits  = ddlraw.fDigits;
70     }
71   return *this;
72 }
73 //____________________________________________________________________________
74
75 AliPMDDDLRawData::~AliPMDDDLRawData()
76 {
77   // Default Destructor
78   //
79
80 }
81
82 //____________________________________________________________________________
83 void AliPMDDDLRawData::WritePMDRawData(TTree *treeD)
84 {
85   // write digits into raw data format
86
87   AliFstream *outfile;
88
89   TBranch *branch = treeD->GetBranch("PMDDigit");
90   if (!branch)
91     {
92       AliError("PMD Digit branch not found");
93       return;
94     }
95   branch->SetAddress(&fDigits);  
96   
97   Int_t   nmodules = (Int_t) treeD->GetEntries();
98   AliDebug(1,Form("Number of modules inside treeD = %d",nmodules));
99
100   const Int_t kDDL          = AliDAQ::NumberOfDdls("PMD");
101
102
103   AliRawDataHeaderSim header;
104   UInt_t sizeRawData = 0;
105   
106   const Int_t kbusSize = 51;
107   const Int_t kSize = 1536;
108   UInt_t buffer[kSize];
109
110   UInt_t busPatch[kbusSize][1536];
111
112   Int_t contentsBus[kbusSize];
113
114   Char_t filename[80];
115
116   Int_t modulePerDDL        = 0;
117   Int_t mmodule = 0;
118   Int_t ddlno;
119   Int_t modulenoddl[12] = {-1,-1,-1,-1,-1,-1,-1,-1,-1,-1,-1,-1};
120
121   for(Int_t iddl = 0; iddl < kDDL; iddl++)
122     {
123       ddlno = iddl;
124       modulePerDDL = fDdlinfo->GetNoOfModulePerDdl(iddl);
125       if (modulePerDDL == 0) continue;
126       for (Int_t im = 0; im < 12; im++)
127         {
128           modulenoddl[im] = fDdlinfo->GetModulesPerDdl(iddl,im);;
129         }
130
131       strncpy(filename,AliDAQ::DdlFileName("PMD",iddl),79);
132       
133       outfile = new AliFstream(filename);
134       
135       // Write the Dummy Data Header into the file
136       Int_t bHPosition = outfile->Tellp();
137       outfile->WriteBuffer((char*)(&header),sizeof(header));
138
139       for (Int_t ibus = 0; ibus < kbusSize; ibus++)
140         {
141           contentsBus[ibus] = 0;
142           for (Int_t ich = 0; ich < kSize; ich++)
143             {
144               busPatch[ibus][ich] = 0;
145             }
146         }
147
148       for(Int_t ium = 0; ium < 12; ium++)
149         {
150           // Extract energy deposition per cell and pack it
151           // in a 32-bit word and returns all the total words
152           // per one unit-module
153           
154           mmodule = modulenoddl[ium];
155           if(mmodule == -1) continue;
156           GetUMDigitsData(treeD, mmodule, iddl, contentsBus, busPatch);
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       dspBlockARDL += 50;
209       dspBlockBRDL += 50;
210
211       // Start writing the DDL file
212
213       AliPMDBlockHeader blHeader;
214       AliPMDDspHeader   dspHeader;
215       AliPMDPatchBusHeader pbusHeader;
216
217       const Int_t kblHLen   = blHeader.GetHeaderLength();
218       const Int_t kdspHLen  = dspHeader.GetHeaderLength();
219       const Int_t kpbusHLen = pbusHeader.GetHeaderLength();
220
221       UInt_t dspRDL = 0;
222       UInt_t dspBlockHeaderWord[8];
223       UInt_t dspHeaderWord[10];
224       UInt_t patchBusHeaderWord[4];
225       Int_t  iskip[5];
226       UInt_t ddlEndWord[2] = {0xDEADFACE, 0xDEADFACE};
227
228       Int_t bknJunk = 0;
229
230
231       for (Int_t iblock = 0; iblock < 2; iblock++)
232         {
233           // DSP Block Header
234           
235           for (Int_t i=0; i<kblHLen; i++)
236             {
237               dspBlockHeaderWord[i] = 0;
238             }
239           if (iblock == 0)
240             {
241               dspBlockHeaderWord[1] = (UInt_t) (dspBlockARDL + kblHLen);
242               dspBlockHeaderWord[2] = (UInt_t) dspBlockARDL;
243             }
244           else if (iblock == 1)
245             {
246               dspBlockHeaderWord[1] = (UInt_t) (dspBlockBRDL + kblHLen);
247               dspBlockHeaderWord[2] = (UInt_t) dspBlockBRDL;
248             }
249
250           outfile->WriteBuffer((char*)dspBlockHeaderWord,kblHLen*sizeof(UInt_t));
251
252           if (iblock == 0)
253             {
254               iskip[0] = 0;
255               iskip[1] = 10;
256               iskip[2] = 20;
257               iskip[3] = 30;
258               iskip[4] = 40;
259             }
260           else if (iblock == 1)
261             {
262               iskip[0] = 5;
263               iskip[1] = 15;
264               iskip[2] = 25;
265               iskip[3] = 35;
266               iskip[4] = 45;
267             }
268
269           for (Int_t idsp = 0; idsp < 5; idsp++)
270             {
271               // DSP Header
272               Int_t dspno = 0;
273               if (iblock == 0)
274                 {
275                   dspno = 2*idsp;
276                   dspRDL = (UInt_t) dsp[dspno];
277                 }
278               else if (iblock == 1)
279                 {
280                   dspno = 2*idsp + 1;
281                   dspRDL = (UInt_t) dsp[dspno];
282                 }
283
284               for (Int_t i=0; i<kdspHLen; i++)
285                 {
286                   dspHeaderWord[i] = 0;
287                 }
288               remainder = dspRDL%2;
289               if (remainder == 1) dspRDL++;
290
291               dspHeaderWord[1] = dspRDL + kdspHLen;
292               dspHeaderWord[2] = dspRDL;
293               dspHeaderWord[3] = dspno;
294               if (remainder == 1) dspHeaderWord[8] = 1; // setting the padding word
295
296
297               outfile->WriteBuffer((char*)dspHeaderWord,kdspHLen*sizeof(UInt_t));
298
299               for (Int_t ibus = 0; ibus < 5; ibus++)
300                 {
301                   // Patch Bus Header
302
303                   Int_t busno = iskip[idsp] + ibus + 1;
304                   Int_t patchbusRDL = contentsBus[busno];
305
306                   if (patchbusRDL > 0)
307                     {
308                       patchBusHeaderWord[0] = 0;
309                       patchBusHeaderWord[1] = (UInt_t) (patchbusRDL + kpbusHLen);
310                       patchBusHeaderWord[2] = (UInt_t) patchbusRDL;
311                       patchBusHeaderWord[3] = (UInt_t) busno;
312                     }
313                   else if (patchbusRDL == 0)
314                     {
315                       patchBusHeaderWord[0] = 0;
316                       patchBusHeaderWord[1] = (UInt_t) kpbusHLen;
317                       patchBusHeaderWord[2] = (UInt_t) 0;
318                       patchBusHeaderWord[3] = (UInt_t) busno;
319                     }
320
321
322                   outfile->WriteBuffer((char*)patchBusHeaderWord,4*sizeof(UInt_t));
323
324                   bknJunk += patchbusRDL;
325
326                   for (Int_t iword = 0; iword < patchbusRDL; iword++)
327                     {
328                       buffer[iword] = busPatch[busno][iword];
329                     }
330                   
331                   outfile->WriteBuffer((char*)buffer,patchbusRDL*sizeof(UInt_t));
332
333                 } // End of patch bus loop
334
335
336               // Adding a padding word if the total words odd
337               if (remainder == 1)
338                 {
339                   UInt_t paddingWord = dspHeader.GetDefaultPaddingWord();
340                   outfile->WriteBuffer((char*)(&paddingWord),sizeof(UInt_t));
341                 }
342             }
343         }
344
345       // Write two extra word at the end of each DDL file
346       outfile->WriteBuffer((char*)ddlEndWord,2*sizeof(UInt_t));
347
348       // Write real data header
349       // take the pointer to the beginning of the data header
350       // write the total number of words per ddl and bring the
351       // pointer to the current file position and close it
352       UInt_t cFPosition = outfile->Tellp();
353       sizeRawData = cFPosition - bHPosition - sizeof(header);
354
355       header.fSize = cFPosition - bHPosition;
356       header.SetAttribute(0);  // valid data
357       outfile->Seekp(bHPosition);
358       outfile->WriteBuffer((char*)(&header),sizeof(header));
359       outfile->Seekp(cFPosition);
360
361       delete outfile;
362     } // DDL Loop over
363
364 }
365 //____________________________________________________________________________
366
367 void AliPMDDDLRawData::GetUMDigitsData(TTree *treeD, Int_t imodule,
368                                        Int_t ddlno,  Int_t *contentsBus,
369                                        UInt_t busPatch[][1536])
370 {
371   // Retrieves digits data UnitModule by UnitModule
372
373   const Int_t kMaxBus = 51;
374
375   UInt_t baseword = 0;
376   UInt_t mcmno = 0, chno = 0;
377   UInt_t adc = 0;
378   Int_t  det = 0, smn = 0, irow = 0, icol = 0;
379   Int_t  parity = 0;
380   
381   //  Int_t totPatchBus = 0, bPatchBus = 0, ePatchBus = 0;
382   //  Int_t ibus = 0, totmcm = 0, rows = 0, cols = 0, rowe = 0, cole = 0;
383   //  Int_t moduleno = 0;
384   Int_t busno = 0;
385   Int_t patchBusNo[kMaxBus], mcmperBus[kMaxBus];
386   Int_t startRowBus[kMaxBus], startColBus[kMaxBus];
387   Int_t endRowBus[kMaxBus], endColBus[kMaxBus];
388
389   Int_t beginPatchBus = -1;
390   Int_t endPatchBus   = -1;
391   for(Int_t i = 0; i < kMaxBus; i++)
392     {
393       patchBusNo[i]  = -1;
394       mcmperBus[i]   = -1;
395       startRowBus[i] = -1;
396       startColBus[i] = -1;
397       endRowBus[i]   = -1;
398       endColBus[i]   = -1;
399     }
400
401   // Fetch the DDL mapping info from the mapping database
402
403   DdlMapping(ddlno, imodule, beginPatchBus, endPatchBus,
404              patchBusNo, mcmperBus, startRowBus, endRowBus,
405              startColBus, endColBus);
406
407   // Read if some chains are off from the ddlinfo database
408
409   Int_t srowoff1[2][24], erowoff1[2][24];
410   Int_t scoloff1[2][24], ecoloff1[2][24];
411   Int_t srowoff2[2][24], erowoff2[2][24];
412   Int_t scoloff2[2][24], ecoloff2[2][24];
413
414   for (Int_t idet = 0; idet < 2; idet++)
415     {
416       for (Int_t im = 0; im < 24; im++)
417         {
418           srowoff1[idet][im] = fDdlinfo->GetStartRowA(idet,im);
419           erowoff1[idet][im] = fDdlinfo->GetEndRowA(idet,im);
420           scoloff1[idet][im] = fDdlinfo->GetStartColA(idet,im);
421           ecoloff1[idet][im] = fDdlinfo->GetEndColA(idet,im);
422           srowoff2[idet][im] = fDdlinfo->GetStartRowB(idet,im);
423           erowoff2[idet][im] = fDdlinfo->GetEndRowB(idet,im);
424           scoloff2[idet][im] = fDdlinfo->GetStartColB(idet,im);
425           ecoloff2[idet][im] = fDdlinfo->GetEndColB(idet,im);
426         }
427     }
428
429   treeD->GetEntry(imodule); 
430   Int_t nentries = fDigits->GetLast();
431   Int_t totword = nentries+1;
432
433   AliPMDdigit *pmddigit = 0x0;
434
435   for (Int_t ient = 0; ient < totword; ient++)
436     {
437       pmddigit = (AliPMDdigit*)fDigits->UncheckedAt(ient);
438       
439       det    = pmddigit->GetDetector();
440       smn    = pmddigit->GetSMNumber();
441       irow   = pmddigit->GetRow();
442       icol   = pmddigit->GetColumn();
443       adc    = (UInt_t) pmddigit->GetADC();
444
445       TransformS2H(smn,irow,icol);
446       
447       // remove the non-existence channels
448
449       //printf("%d %d %d %d\n",det,smn,irow,icol);
450       //printf("--- %d   %d   %d   %d\n",srowoff[det][smn],erowoff[det][smn],
451       //     scoloff[det][smn],ecoloff[det][smn]);
452
453       if (irow >= srowoff1[det][smn] && irow <= erowoff1[det][smn])
454         {
455           if (icol >= scoloff1[det][smn] && icol <= ecoloff1[det][smn])
456             {
457               continue;
458             }
459         }
460       if (irow >= srowoff2[det][smn] && irow <= erowoff2[det][smn])
461         {
462           if (icol >= scoloff2[det][smn] && icol <= ecoloff2[det][smn])
463             {
464               continue;
465             }
466         }
467
468
469       GetMCMCh(imodule, irow, icol, beginPatchBus, endPatchBus,
470                mcmperBus, startRowBus, startColBus,
471                endRowBus, endColBus, busno, mcmno, chno);
472
473       baseword = 0;
474       AliBitPacking::PackWord(adc,baseword,0,11);
475       AliBitPacking::PackWord(chno,baseword,12,17);
476       AliBitPacking::PackWord(mcmno,baseword,18,28);
477       AliBitPacking::PackWord(0,baseword,29,30);
478       parity = ComputeParity(baseword);      // generate the parity bit
479       AliBitPacking::PackWord(parity,baseword,31,31);
480
481       Int_t jj = contentsBus[busno];
482       busPatch[busno][jj] = baseword;
483
484       contentsBus[busno]++;
485     }
486
487 }
488
489 //____________________________________________________________________________
490 void AliPMDDDLRawData::TransformS2H(Int_t smn, Int_t &irow, Int_t &icol)
491 {
492   // Does the Software to Hardware coordinate transformation
493   //
494
495   Int_t  irownew = 0;
496   Int_t  icolnew = 0;
497
498   // First in digits we have all dimension 48x96
499   // Transform into the realistic one, i.e, For SM 0&1 96(row)x48(col)
500   // and for SM 2&3 48(row)x96(col)
501   // 
502   if(smn < 12)
503     {
504       irownew = icol;
505       icolnew = irow;
506     }
507   else if( smn >= 12 && smn < 24)
508     {
509       irownew = irow;
510       icolnew = icol;
511     }
512
513   irow = irownew;
514   icol = icolnew;
515
516 }
517
518
519 //____________________________________________________________________________
520
521 void AliPMDDDLRawData::GetMCMCh(Int_t imodule, Int_t row, Int_t col,
522                                 Int_t beginPatchBus, Int_t endPatchBus,
523                                 Int_t *mcmperBus,
524                                 Int_t *startRowBus, Int_t *startColBus,
525                                 Int_t *endRowBus, Int_t *endColBus,
526                                 Int_t & busno, UInt_t &mcmno, UInt_t &chno)
527 {
528   // This converts row col to hardware channel number
529   // This is the final version of mapping supplied by Mriganka
530
531     UInt_t iCh[16][4];
532
533     static const UInt_t kChDdl01[16][4] = { {6, 4, 5, 7},
534                                            {10, 2, 1, 9},
535                                            {12, 0, 3, 11},
536                                            {14, 8, 13, 15},
537                                            {16, 18, 23, 17},
538                                            {20, 28, 31, 19},
539                                            {22, 30, 29, 21},
540                                            {24, 26, 27, 25},
541                                            {38, 36, 37, 39},
542                                            {42, 34, 33, 41},
543                                            {44, 32, 35, 43},
544                                            {46, 40, 45, 47},
545                                            {48, 50, 55, 49},
546                                            {52, 60, 63, 51},
547                                            {54, 62, 61, 53},
548                                            {56, 58, 59, 57} };
549
550
551     static const UInt_t kChDdl23[16][4] = { {57, 59, 58, 56},
552                                             {53, 61, 62, 54},
553                                             {51, 63, 60, 52},
554                                             {49, 55, 50, 48},
555                                             {47, 45, 40, 46},
556                                             {43, 35, 32, 44},
557                                             {41, 33, 34, 42},
558                                             {39, 37, 36, 38},
559                                             {25, 27, 26, 24},
560                                             {21, 29, 30, 22},
561                                             {19, 31, 28, 20},
562                                             {17, 23, 18, 16},
563                                             {15, 13, 8, 14},
564                                             {11, 3, 0, 12},
565                                             {9, 1, 2, 10},
566                                             {7, 5, 4, 6} };
567     
568     
569     static const UInt_t kChDdl41[16][4] = { {56, 58, 59, 57},
570                                            {54, 62, 61, 53},
571                                            {52, 60, 63, 51},
572                                            {48, 50, 55, 49},
573                                            {46, 40, 45, 47},
574                                            {44, 32, 35, 43},
575                                            {42, 34, 33, 41},
576                                            {38, 36, 37, 39},
577                                            {24, 26, 27, 25},
578                                            {22, 30, 29, 21},
579                                            {20, 28, 31, 19},
580                                            {16, 18, 23, 17},
581                                            {14, 8, 13, 15},
582                                            {12, 0, 3, 11},
583                                            {10, 2, 1, 9},
584                                            {6, 4, 5, 7} };
585
586
587     static const UInt_t kChDdl42[16][4] = { {7, 5, 4, 6},
588                                             {9, 1, 2, 10},
589                                             {11, 3, 0, 12},
590                                             {15, 13, 8, 14},
591                                             {17, 23, 18, 16},
592                                             {19, 31, 28, 20},
593                                             {21, 29, 30, 22},
594                                             {25, 27, 26, 24},
595                                             {39, 37, 36, 38},
596                                             {41, 33, 34, 42},
597                                             {43, 35, 32, 44},
598                                             {47, 45, 40, 46},
599                                             {49, 55, 50, 48},
600                                             {51, 63, 60, 52},
601                                             {53, 61, 62, 54},
602                                             {57, 59, 58, 56} };
603
604
605     static const UInt_t kChDdl51[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
624     static const UInt_t kChDdl52[16][4] = { {56, 58, 59, 57},
625                                             {54, 62, 61, 53},
626                                             {52, 60, 63, 51},
627                                             {48, 50, 55, 49},
628                                             {46, 40, 45, 47},
629                                             {44, 32, 35, 43},
630                                             {42, 34, 33, 41},
631                                             {38, 36, 37, 39},
632                                             {24, 26, 27, 25},
633                                             {22, 30, 29, 21},
634                                             {20, 28, 31, 19},
635                                             {16, 18, 23, 17},
636                                             {14, 8, 13, 15},
637                                             {12, 0, 3, 11},
638                                             {10, 2, 1, 9},
639                                             {6, 4, 5, 7} };
640     
641     
642     for (Int_t i = 0; i < 16; i++)
643       {
644         for (Int_t j = 0; j < 4; j++)
645           {
646             
647             if(imodule < 6)                    iCh[i][j] = kChDdl01[i][j];
648             if(imodule >= 6 && imodule <= 11)  iCh[i][j] = kChDdl01[i][j];
649             if(imodule >= 12 && imodule <= 17) iCh[i][j] = kChDdl23[i][j];
650             if(imodule >= 18 && imodule <= 23) iCh[i][j] = kChDdl23[i][j];
651
652             if(imodule >= 24 && imodule <= 29) iCh[i][j] = kChDdl41[i][j];
653             if(imodule >= 42 && imodule <= 47) iCh[i][j] = kChDdl42[i][j];
654             if(imodule >= 36 && imodule <= 41) iCh[i][j] = kChDdl51[i][j];
655             if(imodule >= 30 && imodule <= 35) iCh[i][j] = kChDdl52[i][j];
656             
657           }
658       }
659
660
661   Int_t irownew = row%16;
662   Int_t icolnew = col%4;
663   
664   chno  = iCh[irownew][icolnew];
665   
666   
667   for (Int_t ibus = beginPatchBus; ibus <= endPatchBus; ibus++)
668     {
669       Int_t srow = startRowBus[ibus];
670       Int_t erow = endRowBus[ibus];
671       Int_t scol = startColBus[ibus];
672       Int_t ecol = endColBus[ibus];
673       Int_t tmcm = mcmperBus[ibus];
674
675       if ((row >= srow && row <= erow) && (col >= scol && col <= ecol))
676         {
677           busno = ibus;
678           
679           // Find out the MCM Number
680           //
681
682           if (imodule < 6)                  mcmno = (col-scol)/4 + 1;
683           if (imodule >= 6 && imodule < 12) mcmno = (col-scol)/4 + 1;
684
685           if (imodule >= 12 && imodule < 18)
686             {
687               icolnew = (col - scol)/4;
688               mcmno = tmcm - icolnew;
689             }         
690           if (imodule >= 18 && imodule < 24)
691             {
692               icolnew = (col - scol)/4;
693               mcmno = tmcm - icolnew;
694             }         
695
696           // DDL = 4
697           if (imodule >= 24 && imodule < 30)
698             {
699
700               //if (tmcm == 24)
701               Int_t rowdiff = endRowBus[ibus] - startRowBus[ibus];
702               if(rowdiff > 16)
703                 {
704                   Int_t midrow = srow + 16;
705                   if(row >= srow && row < midrow)
706                     {
707                       mcmno = 12 + (col-scol)/4 + 1;
708                     }
709                   else if(row >= midrow && row <= erow)
710                     
711                     {
712                       mcmno = (col-scol)/4 + 1;
713                     }
714                 }
715               else if (rowdiff < 16)
716                 {
717                   mcmno = (col-scol)/4 + 1;
718                 }
719             
720             }         
721           if (imodule >= 42 && imodule < 48)
722             {
723               Int_t rowdiff = endRowBus[ibus] - startRowBus[ibus];
724               if(rowdiff > 16)
725                 {
726                   Int_t midrow = srow + 16;
727                   if (row >= midrow && row <= erow)
728                     {
729                       mcmno = 12 + (ecol -col)/4 + 1;
730                     }
731                   else if (row >= srow && row < midrow)
732                     {
733                       mcmno = (ecol - col)/4 + 1;
734                     }
735                 }
736               else if (rowdiff < 16)
737                 {
738                   mcmno = (ecol - col)/4 + 1;
739                 }
740             }         
741
742           // DDL = 5
743           if (imodule >= 30 && imodule < 36)
744             {
745               // CPV plane,  SU Mod = 1, 2 : ddl = 5
746               
747               //if(tmcm == 24)
748               Int_t rowdiff = endRowBus[ibus] - startRowBus[ibus];
749               if(rowdiff > 16)
750                 {
751                   Int_t midrow = srow + 16;
752                   if(row >= srow && row < midrow)
753                     {
754                       mcmno = 12 + (col-scol)/4 + 1;
755                     }
756                   else if(row >= midrow && row <= erow)
757                     {
758                       mcmno = (col-scol)/4 + 1;
759                     }
760                 }
761               else if(rowdiff < 16)
762                 {
763                   mcmno = (col-scol)/4 + 1;
764                 }
765               
766             }
767           if (imodule >= 36 && imodule < 42)
768             {
769               Int_t rowdiff = endRowBus[ibus] - startRowBus[ibus];
770               if(rowdiff > 16)
771                 {
772                   Int_t midrow = srow + 16;
773                   if (row >= midrow && row <= erow)
774                     {
775                       mcmno = 12 + (ecol - col)/4 + 1;
776                     }
777                   else if (row >= srow && row < midrow)
778                     {
779                       mcmno = (ecol - col)/4 + 1;
780                     }
781                 }
782               else if (rowdiff < 16)
783                 {
784                   mcmno = (ecol - col)/4 + 1;
785                 }
786             }
787
788         }
789     }
790
791
792 //____________________________________________________________________________
793
794 Int_t AliPMDDDLRawData::ComputeParity(UInt_t baseword)
795 {
796   // Generate the parity bit
797
798   Int_t count = 0;
799   for(Int_t j=0; j<29; j++)
800     {
801       if (baseword & 0x01 ) count++;
802       baseword >>= 1;
803     }
804   Int_t parity = count%2;
805   return parity;
806 }
807 //____________________________________________________________________________
808 void AliPMDDDLRawData::DdlMapping(Int_t iddl, Int_t imodule,
809                                   Int_t &beginPatchBus, Int_t &endPatchBus,
810                                   Int_t patchBusNo[], Int_t mcmperBus[],
811                                   Int_t startRowBus[], Int_t endRowBus[],
812                                   Int_t startColBus[], Int_t endColBus[])
813 {
814   // DDL Mapping fetching from mapping database
815
816   beginPatchBus = fMapData->GetBeginPatchBus(iddl,imodule);
817   endPatchBus   = fMapData->GetEndPatchBus(iddl,imodule);
818
819   for(Int_t ibus = beginPatchBus; ibus < endPatchBus+1; ibus++)
820     {
821       patchBusNo[ibus]   = ibus;
822       mcmperBus[ibus]    = fMapData->GetMcmperBus(iddl,ibus);
823       startRowBus[ibus]  = fMapData->GetStartRowBus(iddl,ibus);
824       startColBus[ibus]  = fMapData->GetStartColBus(iddl,ibus);
825       endRowBus[ibus]    = fMapData->GetEndRowBus(iddl,ibus);
826       endColBus[ibus]    = fMapData->GetEndColBus(iddl,ibus);
827     }
828
829 }
830 //____________________________________________________________________________
831
832 AliPMDddlinfoData* AliPMDDDLRawData::GetDdlinfoData() const
833 {
834   // The run number will be centralized in AliCDBManager,
835   // you don't need to set it here!
836   AliCDBEntry  *entry = AliCDBManager::Instance()->Get("PMD/Calib/Ddlinfo");
837   
838   if(!entry) AliFatal("ddlinfo object retrieval failed!");
839   
840   AliPMDddlinfoData *ddlinfo = 0;
841   if (entry) ddlinfo = (AliPMDddlinfoData*) entry->GetObject();
842   
843   if (!ddlinfo)  AliFatal("No ddl info data from  database !");
844   
845   return ddlinfo;
846 }
847 //____________________________________________________________________________
848 AliPMDMappingData* AliPMDDDLRawData::GetMappingData() const
849 {
850   // The run number will be centralized in AliCDBManager,
851   // you don't need to set it here!
852   AliCDBEntry  *entry = AliCDBManager::Instance()->Get("PMD/Calib/Mapping");
853   
854   if(!entry) AliFatal("Mapping object retrieval failed!");
855   
856   AliPMDMappingData *mapda = 0;
857   if (entry) mapda = (AliPMDMappingData*) entry->GetObject();
858   
859   if (!mapda)  AliFatal("No mapping data from  database !");
860   
861   return mapda;
862 }
863 //____________________________________________________________________________