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