1 /***************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 #include <Riostream.h>
17 #include <TClonesArray.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 "AliPMDDDLRawData.h"
34 #include "AliFstream.h"
36 ClassImp(AliPMDDDLRawData)
38 AliPMDDDLRawData::AliPMDDDLRawData():
39 fDigits(new TClonesArray("AliPMDdigit", 1000))
41 // Default Constructor
45 //____________________________________________________________________________
46 AliPMDDDLRawData::AliPMDDDLRawData(const AliPMDDDLRawData& ddlraw):
48 fDigits(ddlraw.fDigits)
52 //____________________________________________________________________________
53 AliPMDDDLRawData & AliPMDDDLRawData::operator=(const AliPMDDDLRawData& ddlraw)
58 fDigits = ddlraw.fDigits;
62 //____________________________________________________________________________
64 AliPMDDDLRawData::~AliPMDDDLRawData()
71 //____________________________________________________________________________
72 void AliPMDDDLRawData::WritePMDRawData(TTree *treeD)
74 // write digits into raw data format
78 TBranch *branch = treeD->GetBranch("PMDDigit");
81 AliError("PMD Digit branch not found");
84 branch->SetAddress(&fDigits);
86 Int_t nmodules = (Int_t) treeD->GetEntries();
87 AliDebug(1,Form("Number of modules inside treeD = %d",nmodules));
89 const Int_t kDDL = AliDAQ::NumberOfDdls("PMD");
90 Int_t modulePerDDL = 0;
93 AliRawDataHeaderSim header;
94 UInt_t sizeRawData = 0;
96 const Int_t kbusSize = 51;
97 const Int_t kSize = 1536;
100 UInt_t busPatch[kbusSize][1536];
102 Int_t contentsBus[kbusSize];
108 for(Int_t iddl = 0; iddl < kDDL; iddl++)
110 strcpy(filename,AliDAQ::DdlFileName("PMD",iddl));
112 outfile = new AliFstream(filename);
132 // Write the Dummy Data Header into the file
133 Int_t bHPosition = outfile->Tellp();
134 outfile->WriteBuffer((char*)(&header),sizeof(header));
136 for (Int_t ibus = 0; ibus < kbusSize; ibus++)
138 contentsBus[ibus] = 0;
139 for (Int_t ich = 0; ich < kSize; ich++)
141 busPatch[ibus][ich] = 0;
145 for(Int_t ium = 0; ium < modulePerDDL; ium++)
147 if (iddl == 4 && ium == 6) mmodule = 42;
149 // Extract energy deposition per cell and pack it
150 // in a 32-bit word and returns all the total words
151 // per one unit-module
153 GetUMDigitsData(treeD, mmodule, iddl, contentsBus, busPatch);
162 for (Int_t i = 0; i < 10; i++)
166 for (Int_t ibus=0; ibus < 5; ibus++)
169 if (contentsBus[ij] > 0)
171 dsp[i] += contentsBus[ij];
176 // Add the patch Bus header to the DSP contents
177 dsp[i] += 4*dspBus[i];
180 Int_t dspBlockARDL = 0;
181 Int_t dspBlockBRDL = 0;
183 for (Int_t i = 0; i < 5; i++)
186 Int_t iodd = 2*i + 1;
189 dspBlockARDL += dsp[ieven];
190 Int_t remainder = dsp[ieven]%2;
198 dspBlockBRDL += dsp[iodd];
199 Int_t remainder = dsp[ieven]%2;
207 // Start writing the DDL file
209 AliPMDBlockHeader blHeader;
210 AliPMDDspHeader dspHeader;
211 AliPMDPatchBusHeader pbusHeader;
213 const Int_t kblHLen = blHeader.GetHeaderLength();
214 const Int_t kdspHLen = dspHeader.GetHeaderLength();
215 const Int_t kpbusHLen = pbusHeader.GetHeaderLength();
218 UInt_t dspBlockHeaderWord[8];
219 UInt_t dspHeaderWord[10];
220 UInt_t patchBusHeaderWord[4];
224 for (Int_t iblock = 0; iblock < 2; iblock++)
228 for (Int_t i=0; i<kblHLen; i++)
230 dspBlockHeaderWord[i] = 0;
234 dspBlockHeaderWord[1] = (UInt_t) (dspBlockARDL + kblHLen);
235 dspBlockHeaderWord[2] = (UInt_t) dspBlockARDL;
237 else if (iblock == 1)
239 dspBlockHeaderWord[1] = (UInt_t) (dspBlockBRDL + kblHLen);
240 dspBlockHeaderWord[2] = (UInt_t) dspBlockBRDL;
243 outfile->WriteBuffer((char*)dspBlockHeaderWord,kblHLen*sizeof(UInt_t));
253 else if (iblock == 1)
262 for (Int_t idsp = 0; idsp < 5; idsp++)
269 dspRDL = (UInt_t) dsp[dspno];
271 else if (iblock == 1)
274 dspRDL = (UInt_t) dsp[dspno];
277 for (Int_t i=0; i<kdspHLen; i++)
279 dspHeaderWord[i] = 0;
281 Int_t remainder = dspRDL%2;
282 if (remainder == 1) dspRDL++;
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->WriteBuffer((char*)dspHeaderWord,kdspHLen*sizeof(UInt_t));
290 for (Int_t ibus = 0; ibus < 5; ibus++)
293 // BKN - just added 1
294 Int_t busno = iskip[idsp] + ibus + 1;
295 Int_t patchbusRDL = contentsBus[busno];
299 patchBusHeaderWord[0] = 0;
300 patchBusHeaderWord[1] = (UInt_t) (patchbusRDL + kpbusHLen);
301 patchBusHeaderWord[2] = (UInt_t) patchbusRDL;
302 patchBusHeaderWord[3] = (UInt_t) busno;
304 else if (patchbusRDL == 0)
306 patchBusHeaderWord[0] = 0;
307 patchBusHeaderWord[1] = (UInt_t) kpbusHLen;
308 patchBusHeaderWord[2] = (UInt_t) 0;
309 patchBusHeaderWord[3] = (UInt_t) busno;
313 outfile->WriteBuffer((char*)patchBusHeaderWord,4*sizeof(UInt_t));
316 for (Int_t iword = 0; iword < patchbusRDL; iword++)
318 buffer[iword] = busPatch[busno][iword];
321 outfile->WriteBuffer((char*)buffer,patchbusRDL*sizeof(UInt_t));
323 } // End of patch bus loop
326 // Adding a padding word if the total words odd
329 UInt_t paddingWord = dspHeader.GetDefaultPaddingWord();
330 outfile->WriteBuffer((char*)(&paddingWord),sizeof(UInt_t));
335 // Write real data header
336 // take the pointer to the beginning of the data header
337 // write the total number of words per ddl and bring the
338 // pointer to the current file position and close it
339 UInt_t cFPosition = outfile->Tellp();
340 sizeRawData = cFPosition - bHPosition - sizeof(header);
342 header.fSize = cFPosition - bHPosition;
343 header.SetAttribute(0); // valid data
344 outfile->Seekp(bHPosition);
345 outfile->WriteBuffer((char*)(&header),sizeof(header));
346 outfile->Seekp(cFPosition);
353 //____________________________________________________________________________
354 void AliPMDDDLRawData::GetUMDigitsData(TTree *treeD, Int_t imodule,
355 Int_t ddlno, Int_t *contentsBus,
356 UInt_t busPatch[][1536])
358 // Retrives digits data UnitModule by UnitModule
363 Int_t det, smn, irow, icol;
366 const Int_t kMaxBus = 51; // BKN
367 Int_t totPatchBus, bPatchBus, ePatchBus;
368 Int_t ibus, totmcm, rows, cols, rowe, cole;
371 Int_t patchBusNo[kMaxBus], mcmperBus[kMaxBus];
372 Int_t startRowBus[kMaxBus], startColBus[kMaxBus];
373 Int_t endRowBus[kMaxBus], endColBus[kMaxBus];
375 Int_t beginPatchBus = -1;
376 Int_t endPatchBus = -1;
377 for(Int_t i = 0; i < kMaxBus; i++)
386 Int_t modulePerDDL = 0;
391 else if (ddlno == 4 || ddlno == 5)
396 TString fileName(gSystem->Getenv("ALICE_ROOT"));
400 fileName += "/PMD/PMD_Mapping_ddl0.dat";
404 fileName += "/PMD/PMD_Mapping_ddl1.dat";
408 fileName += "/PMD/PMD_Mapping_ddl2.dat";
412 fileName += "/PMD/PMD_Mapping_ddl3.dat";
416 fileName += "/PMD/PMD_Mapping_ddl4.dat";
420 fileName += "/PMD/PMD_Mapping_ddl5.dat";
424 infile.open(fileName.Data(), ios::in); // ascii file
426 AliError(Form("Could not read the mapping file for DDL No = %d",ddlno));
428 for (Int_t im = 0; im < modulePerDDL; im++)
431 infile >> totPatchBus >> bPatchBus >> ePatchBus;
433 if (moduleno == imodule)
435 beginPatchBus = bPatchBus;
436 endPatchBus = ePatchBus;
439 for(Int_t i=0; i<totPatchBus; i++)
441 infile >> ibus >> totmcm >> rows >> rowe >> cols >> cole;
443 if (moduleno == imodule)
445 patchBusNo[ibus] = ibus;
446 mcmperBus[ibus] = totmcm;
447 startRowBus[ibus] = rows;
448 startColBus[ibus] = cols;
449 endRowBus[ibus] = rowe;
450 endColBus[ibus] = cole;
458 treeD->GetEntry(imodule);
459 Int_t nentries = fDigits->GetLast();
460 Int_t totword = nentries+1;
462 AliPMDdigit *pmddigit;
464 for (Int_t ient = 0; ient < totword; ient++)
466 pmddigit = (AliPMDdigit*)fDigits->UncheckedAt(ient);
468 det = pmddigit->GetDetector();
469 smn = pmddigit->GetSMNumber();
470 irow = pmddigit->GetRow();
471 icol = pmddigit->GetColumn();
472 adc = (UInt_t) pmddigit->GetADC();
474 TransformS2H(smn,irow,icol);
476 GetMCMCh(ddlno, smn, irow, icol, beginPatchBus, endPatchBus,
477 mcmperBus, startRowBus, startColBus,
478 endRowBus, endColBus, busno, mcmno, chno);
481 AliBitPacking::PackWord(adc,baseword,0,11);
482 AliBitPacking::PackWord(chno,baseword,12,17);
483 AliBitPacking::PackWord(mcmno,baseword,18,28);
484 AliBitPacking::PackWord(0,baseword,29,30);
485 parity = ComputeParity(baseword); // generate the parity bit
486 AliBitPacking::PackWord(parity,baseword,31,31);
488 Int_t jj = contentsBus[busno];
489 busPatch[busno][jj] = baseword;
491 contentsBus[busno]++;
496 //____________________________________________________________________________
497 void AliPMDDDLRawData::TransformS2H(Int_t smn, Int_t &irow, Int_t &icol)
499 // Does the Software to Hardware coordinate transformation
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)
514 else if( smn >= 12 && smn < 24)
526 //____________________________________________________________________________
528 void AliPMDDDLRawData::GetMCMCh(Int_t ddlno, Int_t smn, Int_t row, Int_t col,
529 Int_t beginPatchBus, Int_t endPatchBus,
531 Int_t *startRowBus, Int_t *startColBus,
532 Int_t *endRowBus, Int_t *endColBus,
533 Int_t & busno, UInt_t &mcmno, UInt_t &chno)
535 // This converts row col to hardware channel number
536 // This is the final version of mapping supplied by Mriganka
540 // ChMap(ddlno, smn, iCh);
543 static const UInt_t kChDdl01[16][4] = { {6, 4, 5, 7},
561 static const UInt_t kChDdl23[16][4] = { {57, 59, 58, 56},
579 static const UInt_t kChDdl41[16][4] = { {56, 58, 59, 57},
597 static const UInt_t kChDdl42[16][4] = { {7, 5, 4, 6},
615 static const UInt_t kChDdl51[16][4] = { {7, 5, 4, 6},
634 static const UInt_t kChDdl52[16][4] = { {56, 58, 59, 57},
652 for (Int_t i = 0; i < 16; i++)
654 for (Int_t j = 0; j < 4; j++)
656 if (ddlno == 0 || ddlno == 1) iCh[i][j] = kChDdl01[i][j];
657 if (ddlno == 2 || ddlno == 3) iCh[i][j] = kChDdl23[i][j];
659 if (ddlno == 4 && smn < 6) iCh[i][j] = kChDdl41[i][j];
660 if (ddlno == 4 && (smn >= 18 && smn < 24))iCh[i][j] = kChDdl42[i][j];
661 if (ddlno == 5 && (smn >= 12 && smn < 18))iCh[i][j] = kChDdl51[i][j];
662 if (ddlno == 5 && (smn >= 6 && smn < 12))iCh[i][j] = kChDdl52[i][j];
667 Int_t irownew = row%16;
668 Int_t icolnew = col%4;
670 chno = iCh[irownew][icolnew];
673 for (Int_t ibus = beginPatchBus; ibus <= endPatchBus; ibus++)
675 Int_t srow = startRowBus[ibus];
676 Int_t erow = endRowBus[ibus];
677 Int_t scol = startColBus[ibus];
678 Int_t ecol = endColBus[ibus];
679 Int_t tmcm = mcmperBus[ibus];
680 if ((row >= srow && row <= erow) && (col >= scol && col <= ecol))
684 // Find out the MCM Number
687 if (ddlno == 0 || ddlno == 1)
689 // PRE plane, SU Mod = 0, 1
690 mcmno = (col-scol)/4 + 1;
693 else if (ddlno == 2 || ddlno == 3)
695 // PRE plane, SU Mod = 2, 3
696 Int_t icolnew = (col - scol)/4;
697 mcmno = tmcm - icolnew;
699 else if (ddlno == 4 )
701 // CPV plane, SU Mod = 0, 3 : ddl = 4
705 Int_t midrow = srow + 16;
706 if(row >= srow && row < midrow)
708 mcmno = 12 + (col-scol)/4 + 1;
710 else if(row >= midrow && row <= erow)
713 mcmno = (col-scol)/4 + 1;
719 Int_t rowdiff = endRowBus[ibus] - startRowBus[ibus];
722 Int_t midrow = srow + 16;
723 if (row >= midrow && row <= erow)
725 mcmno = 12 + (ecol -col)/4 + 1;
727 else if (row >= srow && row < midrow)
729 mcmno = (ecol - col)/4 + 1;
732 else if (rowdiff < 16)
734 mcmno = (ecol - col)/4 + 1;
738 else if ( ddlno == 5)
740 // CPV plane, SU Mod = 1, 2 : ddl = 5
744 Int_t midrow = srow + 16;
745 if(row >= srow && row < midrow)
747 mcmno = 12 + (col-scol)/4 + 1;
749 else if(row >= midrow && row <= erow)
751 mcmno = (col-scol)/4 + 1;
757 Int_t rowdiff = endRowBus[ibus] - startRowBus[ibus];
760 Int_t midrow = srow + 16;
761 if (row >= midrow && row <= erow)
763 mcmno = 12 + (ecol - col)/4 + 1;
765 else if (row >= srow && row < midrow)
767 mcmno = (ecol - col)/4 + 1;
770 else if (rowdiff < 16)
772 mcmno = (ecol - col)/4 + 1;
780 //____________________________________________________________________________
783 void ChMap(Int_t ddlno, Int_t smn, UInt_t iCh[][4])
785 // Row and Col converted to Channel
787 static const UInt_t kChDdl01[16][4] = { {6, 4, 5, 7},
805 static const UInt_t kChDdl23[16][4] = { {57, 59, 58, 56},
823 static const UInt_t kChDdl41[16][4] = { {56, 58, 59, 57},
841 static const UInt_t kChDdl42[16][4] = { {7, 5, 4, 6},
859 static const UInt_t kChDdl51[16][4] = { {7, 5, 4, 6},
878 static const UInt_t kChDdl52[16][4] = { {56, 58, 59, 57},
896 for (Int_t i = 0; i < 16; i++)
898 for (Int_t j = 0; j < 4; j++)
900 if (ddlno == 0 || ddlno == 1) iCh[i][j] = kChDdl01[i][j];
901 if (ddlno == 2 || ddlno == 3) iCh[i][j] = kChDdl23[i][j];
903 if (ddlno == 4 && smn < 6) iCh[i][j] = kChDdl41[i][j];
904 if (ddlno == 4 && (smn >= 18 && smn < 24))iCh[i][j] = kChDdl42[i][j];
905 if (ddlno == 5 && (smn >= 12 && smn < 18))iCh[i][j] = kChDdl51[i][j];
906 if (ddlno == 0 && (smn >= 6 && smn < 12))iCh[i][j] = kChDdl52[i][j];
912 //____________________________________________________________________________
916 Int_t AliPMDDDLRawData::ComputeParity(UInt_t baseword)
918 // Generate the parity bit
921 for(Int_t j=0; j<29; j++)
923 if (baseword & 0x01 ) count++;
926 Int_t parity = count%2;
930 //____________________________________________________________________________