]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PMD/AliPMDDDLRawData.cxx
minor bugfix in argument scanning
[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 "AliRawDataHeader.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   AliRawDataHeader header;
94   UInt_t sizeRawData = 0;
95   
96   const Int_t kSize = 1536;
97   UInt_t buffer[kSize];
98
99   UInt_t busPatch[50][1536];
100
101   Int_t contentsBus[50];
102
103   Char_t filename[80];
104
105
106   Int_t mmodule = 0;
107   for(Int_t iddl = 0; iddl < kDDL; iddl++)
108     {
109       strcpy(filename,AliDAQ::DdlFileName("PMD",iddl));
110       
111       outfile = new AliFstream(filename);
112       
113       if (iddl < 4)
114         {
115           modulePerDDL = 6;
116           mmodule = 6*iddl;
117         }
118       else if (iddl == 4)
119         {
120           modulePerDDL = 12;
121           mmodule = 24;
122         }
123       else if (iddl == 5)
124         {
125           modulePerDDL = 12;
126           mmodule = 30;
127         }
128
129
130
131       // Write the Dummy Data Header into the file
132       Int_t bHPosition = outfile->Tellp();
133       outfile->WriteBuffer((char*)(&header),sizeof(header));
134
135       for (Int_t ibus = 0; ibus < 50; ibus++)
136         {
137           contentsBus[ibus] = 0;
138           for (Int_t ich = 0; ich < 1536; ich++)
139             {
140               busPatch[ibus][ich] = 0;
141             }
142         }
143
144       for(Int_t ium = 0; ium < modulePerDDL; ium++)
145         {
146           if (iddl == 4 && ium == 6) mmodule = 42;
147
148           // Extract energy deposition per cell and pack it
149           // in a 32-bit word and returns all the total words
150           // per one unit-module
151           
152           GetUMDigitsData(treeD, mmodule, iddl, contentsBus, busPatch);
153           mmodule++;
154         }
155
156       
157
158       Int_t ij = 0;
159       Int_t dsp[10];
160       Int_t dspBus[10];
161       for (Int_t i = 0; i < 10; i++)
162         {
163           dsp[i] = 0;
164           dspBus[i] = 0;
165           for (Int_t ibus=0; ibus < 5; ibus++)
166             {
167               if (contentsBus[ij] > 0)
168                 {
169                   dsp[i] += contentsBus[ij];
170                   dspBus[i]++;
171                 }
172               ij++;
173             }
174           // Add the patch Bus header to the DSP contents
175           dsp[i] += 4*dspBus[i];
176         }
177
178       Int_t dspBlockARDL    = 0;
179       Int_t dspBlockBRDL    = 0;
180
181       for (Int_t i = 0; i < 5; i++)
182         {
183           Int_t ieven = 2*i;
184           Int_t iodd  = 2*i + 1;
185           if (dsp[ieven] > 0)
186             {
187               dspBlockARDL += dsp[ieven];
188               Int_t remainder = dsp[ieven]%2;
189               if (remainder == 1)
190                 {
191                   dspBlockARDL++;
192                 }
193             }
194           if (dsp[iodd] > 0)
195             {
196               dspBlockBRDL += dsp[iodd];
197               Int_t remainder = dsp[ieven]%2;
198               if (remainder == 1)
199                 {
200                   dspBlockARDL++;
201                 }
202             }
203         }
204       
205       // Start writing the DDL file
206
207       AliPMDBlockHeader blHeader;
208       AliPMDDspHeader   dspHeader;
209       AliPMDPatchBusHeader pbusHeader;
210
211       const Int_t kblHLen   = blHeader.GetHeaderLength();
212       const Int_t kdspHLen  = dspHeader.GetHeaderLength();
213       const Int_t kpbusHLen = pbusHeader.GetHeaderLength();
214
215       UInt_t dspRDL = 0;
216       UInt_t dspBlockHeaderWord[8];
217       UInt_t dspHeaderWord[10];
218       UInt_t patchBusHeaderWord[4];
219       Int_t iskip[5];
220
221
222       for (Int_t iblock = 0; iblock < 2; iblock++)
223         {
224           // DSP Block Header
225           
226           for (Int_t i=0; i<kblHLen; i++)
227             {
228               dspBlockHeaderWord[i] = 0;
229             }
230           if (iblock == 0)
231             {
232               dspBlockHeaderWord[1] = (UInt_t) (dspBlockARDL + kblHLen);
233               dspBlockHeaderWord[2] = (UInt_t) dspBlockARDL;
234             }
235           else if (iblock == 1)
236             {
237               dspBlockHeaderWord[1] = (UInt_t) (dspBlockBRDL + kblHLen);
238               dspBlockHeaderWord[2] = (UInt_t) dspBlockBRDL;
239             }
240
241           outfile->WriteBuffer((char*)dspBlockHeaderWord,kblHLen*sizeof(UInt_t));
242
243           if (iblock == 0)
244             {
245               iskip[0] = 0;
246               iskip[1] = 10;
247               iskip[2] = 20;
248               iskip[3] = 30;
249               iskip[4] = 40;
250             }
251           else if (iblock == 1)
252             {
253               iskip[0] = 5;
254               iskip[1] = 15;
255               iskip[2] = 25;
256               iskip[3] = 35;
257               iskip[4] = 45;
258             }
259
260           for (Int_t idsp = 0; idsp < 5; idsp++)
261             {
262               // DSP Header
263               Int_t dspno = 0;
264               if (iblock == 0)
265                 {
266                   dspno = 2*idsp;
267                   dspRDL = (UInt_t) dsp[dspno];
268                 }
269               else if (iblock == 1)
270                 {
271                   dspno = 2*idsp + 1;
272                   dspRDL = (UInt_t) dsp[dspno];
273                 }
274
275               for (Int_t i=0; i<kdspHLen; i++)
276                 {
277                   dspHeaderWord[i] = 0;
278                 }
279               Int_t remainder = dspRDL%2;
280               if (remainder == 1) dspRDL++;
281
282               dspHeaderWord[1] = dspRDL + kdspHLen;
283               dspHeaderWord[2] = dspRDL;
284               dspHeaderWord[3] = dspno;
285               if (remainder == 1) dspHeaderWord[8] = 1; // setting the padding word
286               outfile->WriteBuffer((char*)dspHeaderWord,kdspHLen*sizeof(UInt_t));
287
288               for (Int_t ibus = 0; ibus < 5; ibus++)
289                 {
290                   // Patch Bus Header
291                   Int_t busno = iskip[idsp] + ibus;
292                   Int_t patchbusRDL = contentsBus[busno];
293
294                   if (patchbusRDL > 0)
295                     {
296                       patchBusHeaderWord[0] = 0;
297                       patchBusHeaderWord[1] = (UInt_t) (patchbusRDL + kpbusHLen);
298                       patchBusHeaderWord[2] = (UInt_t) patchbusRDL;
299                       patchBusHeaderWord[3] = (UInt_t) busno;
300                     }
301                   else if (patchbusRDL == 0)
302                     {
303                       patchBusHeaderWord[0] = 0;
304                       patchBusHeaderWord[1] = (UInt_t) kpbusHLen;
305                       patchBusHeaderWord[2] = (UInt_t) 0;
306                       patchBusHeaderWord[3] = (UInt_t) busno;
307                     }
308
309
310                   outfile->WriteBuffer((char*)patchBusHeaderWord,4*sizeof(UInt_t));
311
312
313                   for (Int_t iword = 0; iword < patchbusRDL; iword++)
314                     {
315                       buffer[iword] = busPatch[busno][iword];
316                     }
317                   
318                   outfile->WriteBuffer((char*)buffer,patchbusRDL*sizeof(UInt_t));
319
320                 } // End of patch bus loop
321
322
323               // Adding a padding word if the total words odd
324               if (remainder == 1)
325                 {
326                   UInt_t paddingWord = dspHeader.GetDefaultPaddingWord();
327                   outfile->WriteBuffer((char*)(&paddingWord),sizeof(UInt_t));
328                 }
329             }
330         }
331
332       // Write real data header
333       // take the pointer to the beginning of the data header
334       // write the total number of words per ddl and bring the
335       // pointer to the current file position and close it
336       UInt_t cFPosition = outfile->Tellp();
337       sizeRawData = cFPosition - bHPosition - sizeof(header);
338
339       header.fSize = cFPosition - bHPosition;
340       header.SetAttribute(0);  // valid data
341       outfile->Seekp(bHPosition);
342       outfile->WriteBuffer((char*)(&header),sizeof(header));
343       outfile->Seekp(cFPosition);
344
345       delete outfile;
346     } // DDL Loop over
347
348
349 }
350 //____________________________________________________________________________
351 void AliPMDDDLRawData::GetUMDigitsData(TTree *treeD, Int_t imodule,
352                                        Int_t ddlno, Int_t *contentsBus,
353                                        UInt_t busPatch[][1536])
354 {
355   // Retrives digits data UnitModule by UnitModule
356
357   UInt_t baseword;
358   UInt_t mcmno, chno;
359   UInt_t adc;
360   Int_t  det, smn, irow, icol;
361   Int_t  parity;
362   
363   const Int_t kMaxBus = 50;
364   Int_t totPatchBus, bPatchBus, ePatchBus;
365   Int_t ibus, totmcm, rows, cols, rowe, cole;
366   Int_t moduleno;
367   Int_t busno = 0;
368   Int_t patchBusNo[kMaxBus], mcmperBus[kMaxBus];
369   Int_t startRowBus[kMaxBus], startColBus[kMaxBus];
370   Int_t endRowBus[kMaxBus], endColBus[kMaxBus];
371
372   Int_t beginPatchBus = -1;
373   Int_t endPatchBus   = -1;
374   for(Int_t i = 0; i < kMaxBus; i++)
375     {
376       patchBusNo[i]  = -1;
377       mcmperBus[i]   = -1;
378       startRowBus[i] = -1;
379       startColBus[i] = -1;
380       endRowBus[i]   = -1;
381       endColBus[i]   = -1;
382     }
383   Int_t modulePerDDL = 0;
384   if (ddlno < 4)
385     {
386       modulePerDDL = 6;
387     }
388   else if (ddlno == 4 || ddlno == 5)
389     {
390       modulePerDDL = 12;
391     }
392
393   TString fileName(gSystem->Getenv("ALICE_ROOT"));
394
395   if(ddlno == 0)
396     {
397       fileName += "/PMD/PMD_Mapping_ddl0.dat";
398     }
399   else if(ddlno == 1)
400     {
401       fileName += "/PMD/PMD_Mapping_ddl1.dat";
402     }
403   else if(ddlno == 2)
404     {
405       fileName += "/PMD/PMD_Mapping_ddl2.dat";
406     }
407   else if(ddlno == 3)
408     {
409       fileName += "/PMD/PMD_Mapping_ddl3.dat";
410     }
411   else if(ddlno == 4)
412     {
413       fileName += "/PMD/PMD_Mapping_ddl4.dat";
414     }
415   else if(ddlno == 5)
416     {
417       fileName += "/PMD/PMD_Mapping_ddl5.dat";
418     }
419
420   ifstream infile;
421   infile.open(fileName.Data(), ios::in); // ascii file
422   if(!infile)
423     AliError(Form("Could not read the mapping file for DDL No = %d",ddlno));
424
425   for (Int_t im = 0; im < modulePerDDL; im++)
426     {
427       infile >> moduleno;
428       infile >> totPatchBus >> bPatchBus >> ePatchBus;
429
430       if (moduleno == imodule)
431         {
432           beginPatchBus = bPatchBus;
433           endPatchBus   = ePatchBus;
434         }
435       
436       for(Int_t i=0; i<totPatchBus; i++)
437         {
438           infile >> ibus >> totmcm >> rows >> rowe >> cols >> cole;
439
440           if (moduleno == imodule)
441             {
442               patchBusNo[ibus]   = ibus;
443               mcmperBus[ibus]    = totmcm;
444               startRowBus[ibus]  = rows;
445               startColBus[ibus]  = cols;
446               endRowBus[ibus]    = rowe;
447               endColBus[ibus]    = cole;
448             }
449         }
450
451     }
452
453   infile.close();
454
455   treeD->GetEntry(imodule); 
456   Int_t nentries = fDigits->GetLast();
457   Int_t totword = nentries+1;
458
459   AliPMDdigit *pmddigit;
460
461   for (Int_t ient = 0; ient < totword; ient++)
462     {
463       pmddigit = (AliPMDdigit*)fDigits->UncheckedAt(ient);
464       
465       det    = pmddigit->GetDetector();
466       smn    = pmddigit->GetSMNumber();
467       irow   = pmddigit->GetRow();
468       icol   = pmddigit->GetColumn();
469       adc    = (UInt_t) pmddigit->GetADC();
470
471       TransformS2H(smn,irow,icol);
472
473       GetMCMCh(ddlno, irow, icol, beginPatchBus, endPatchBus,
474                mcmperBus, startRowBus, startColBus,
475                endRowBus, endColBus, busno, mcmno, chno);
476
477       baseword = 0;
478       AliBitPacking::PackWord(adc,baseword,0,11);
479       AliBitPacking::PackWord(chno,baseword,12,17);
480       AliBitPacking::PackWord(mcmno,baseword,18,28);
481       AliBitPacking::PackWord(0,baseword,29,30);
482       parity = ComputeParity(baseword);      // generate the parity bit
483       AliBitPacking::PackWord(parity,baseword,31,31);
484
485       Int_t jj = contentsBus[busno];
486       busPatch[busno][jj] = baseword;
487
488       contentsBus[busno]++;
489     }
490
491 }
492
493 //____________________________________________________________________________
494 void AliPMDDDLRawData::TransformS2H(Int_t smn, Int_t &irow, Int_t &icol)
495 {
496   // Does the Software to Hardware coordinate transformation
497   //
498
499   Int_t  irownew = 0;
500   Int_t  icolnew = 0;
501
502   // First in digits we have all dimension 48x96
503   // Transform into the realistic one, i.e, For SM 0&1 96(row)x48(col)
504   // and for SM 2&3 48(row)x96(col)
505   // 
506   if(smn < 12)
507     {
508       irownew = icol;
509       icolnew = irow;
510     }
511   else if( smn >= 12 && smn < 24)
512     {
513       irownew = irow;
514       icolnew = icol;
515     }
516
517   irow = irownew;
518   icol = icolnew;
519
520 }
521
522
523 //____________________________________________________________________________
524
525 void AliPMDDDLRawData::GetMCMCh(Int_t ddlno, Int_t row, Int_t col,
526                                 Int_t beginPatchBus, Int_t endPatchBus,
527                                 Int_t *mcmperBus,
528                                 Int_t *startRowBus, Int_t *startColBus,
529                                 Int_t *endRowBus, Int_t *endColBus,
530                                 Int_t & busno, UInt_t &mcmno, UInt_t &chno)
531 {
532   // This converts row col to hardware channel number
533   // This is the final version of mapping supplied by Mriganka
534
535   static const UInt_t kCh[16][4] = { {56, 58, 59, 57},
536                                      {54, 62, 61, 53},
537                                      {52, 60, 63, 51},
538                                      {48, 50, 55, 49},
539                                      {46, 40, 45, 47},
540                                      {44, 32, 35, 43},
541                                      {42, 34, 33, 41},
542                                      {38, 36, 37, 39},
543                                      {24, 26, 27, 25},
544                                      {22, 30, 29, 21},
545                                      {20, 28, 31, 19},
546                                      {16, 18, 23, 17},
547                                      {14, 8, 13, 15},
548                                      {12, 0, 3, 11},
549                                      {10, 2, 1, 9},
550                                      {6, 4, 5, 7} };
551   
552   Int_t irownew = row%16;
553   Int_t icolnew = col%4;
554   
555   chno  = kCh[irownew][icolnew];
556   
557   
558   for (Int_t ibus = beginPatchBus; ibus <= endPatchBus; ibus++)
559     {
560       Int_t srow = startRowBus[ibus];
561       Int_t erow = endRowBus[ibus];
562       Int_t scol = startColBus[ibus];
563       Int_t ecol = endColBus[ibus];
564       Int_t tmcm = mcmperBus[ibus];
565       if ((row >= srow && row <= erow) && (col >= scol && col <= ecol))
566         {
567           busno = ibus;
568           
569           // Find out the MCM Number
570           //
571           
572           if (ddlno == 0 || ddlno == 1)
573             {
574               // PRE plane, SU Mod = 0, 1
575               mcmno = (col-scol)/4;
576               
577             }
578           else if (ddlno == 2 || ddlno == 3)
579             {
580               // PRE plane,  SU Mod = 2, 3
581               Int_t icolnew = (col - scol)/4;
582               mcmno = tmcm - 1 - icolnew;
583             }
584           else if (ddlno == 4 )
585             {
586               // CPV plane,  SU Mod = 0, 3 : ddl = 4
587               
588               if(ibus <= 17)
589                 {
590                   Int_t midrow = srow + 16;
591                   if(row >= srow && row < midrow)
592                     {
593                       mcmno = 12 + (col-scol)/4;
594                     }
595                   else if(row >= midrow && row <= erow)
596                   
597                     {
598                       mcmno = (col-scol)/4;
599                     }
600                 }
601               
602               else if (ibus > 17)
603                 {
604                   Int_t rowdiff = endRowBus[ibus] - startRowBus[ibus];
605                   if(rowdiff > 16)
606                     {
607                       Int_t midrow = srow + 16;
608                       if (row >= midrow && row <= erow)
609                         {
610                           mcmno = 12 + (ecol -col)/4;
611                         }
612                       else if (row >= srow && row < midrow)
613                         {
614                           mcmno = (ecol - col)/4;
615                         }
616                     }
617                   else if (rowdiff < 16) 
618                     {
619                       mcmno = (ecol - col)/4;
620                     } 
621                 }
622             }
623           else if ( ddlno == 5)
624             {
625               // CPV plane,  SU Mod = 0, 3 : ddl = 4
626               
627               if(ibus <= 17)
628                 {
629                   Int_t midrow = srow + 16;
630                   if(row >= srow && row < midrow)
631                     {
632                       mcmno = 12 + (col-scol)/4;
633                     }
634                   else if(row >= midrow && row <= erow)
635                     {
636                       mcmno = (col-scol)/4;
637                     }
638                 }
639               
640               else if (ibus > 17)
641                 {
642                   Int_t rowdiff = endRowBus[ibus] - startRowBus[ibus];
643                   if(rowdiff > 16)
644                     {
645                       Int_t midrow = srow + 16;
646                       if (row >= midrow && row <= erow)
647                         {
648                           mcmno = 12 + (ecol -col)/4;
649                         }
650                       else if (row >= srow && row < midrow)
651                         {
652                           mcmno = (ecol - col)/4;
653                         }
654                     }
655                   else if (rowdiff < 16) 
656                     {
657                       mcmno = (ecol - col)/4;
658                     } 
659                 }
660             }
661         }
662     }
663
664
665 //____________________________________________________________________________
666
667 Int_t AliPMDDDLRawData::ComputeParity(UInt_t baseword)
668 {
669   // Generate the parity bit
670
671   Int_t count = 0;
672   for(Int_t j=0; j<29; j++)
673     {
674       if (baseword & 0x01 ) count++;
675       baseword >>= 1;
676     }
677   Int_t parity = count%2;
678   return parity;
679 }
680
681 //____________________________________________________________________________