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