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