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 "AliPMDddlinfoData.h"
33 #include "AliPMDMappingData.h"
34 #include "AliPMDDDLRawData.h"
36 #include "AliFstream.h"
37 #include "AliCDBManager.h"
38 #include "AliCDBStorage.h"
39 #include "AliCDBEntry.h"
41 ClassImp(AliPMDDDLRawData)
43 AliPMDDDLRawData::AliPMDDDLRawData():
44 fDdlinfo(GetDdlinfoData()),
45 fMapData(GetMappingData()),
46 fDigits(new TClonesArray("AliPMDdigit", 1000))
48 // Default Constructor
52 //____________________________________________________________________________
53 AliPMDDDLRawData::AliPMDDDLRawData(const AliPMDDDLRawData& ddlraw):
55 fDdlinfo(ddlraw.fDdlinfo),
56 fMapData(ddlraw.fMapData),
57 fDigits(ddlraw.fDigits)
61 //____________________________________________________________________________
62 AliPMDDDLRawData & AliPMDDDLRawData::operator=(const AliPMDDDLRawData& ddlraw)
67 fDdlinfo = ddlraw.fDdlinfo;
68 fMapData = ddlraw.fMapData;
69 fDigits = ddlraw.fDigits;
73 //____________________________________________________________________________
75 AliPMDDDLRawData::~AliPMDDDLRawData()
82 //____________________________________________________________________________
83 void AliPMDDDLRawData::WritePMDRawData(TTree *treeD)
85 // write digits into raw data format
89 TBranch *branch = treeD->GetBranch("PMDDigit");
92 AliError("PMD Digit branch not found");
95 branch->SetAddress(&fDigits);
97 Int_t nmodules = (Int_t) treeD->GetEntries();
98 AliDebug(1,Form("Number of modules inside treeD = %d",nmodules));
100 const Int_t kDDL = AliDAQ::NumberOfDdls("PMD");
103 AliRawDataHeaderSim header;
104 UInt_t sizeRawData = 0;
106 const Int_t kbusSize = 51;
107 const Int_t kSize = 1536;
108 UInt_t buffer[kSize];
110 UInt_t busPatch[kbusSize][1536];
112 Int_t contentsBus[kbusSize];
116 Int_t modulePerDDL = 0;
119 Int_t modulenoddl[12] = {-1,-1,-1,-1,-1,-1,-1,-1,-1,-1,-1,-1};
121 for(Int_t iddl = 0; iddl < kDDL; iddl++)
124 modulePerDDL = fDdlinfo->GetNoOfModulePerDdl(iddl);
125 if (modulePerDDL == 0) continue;
126 for (Int_t im = 0; im < 12; im++)
128 modulenoddl[im] = fDdlinfo->GetModulesPerDdl(iddl,im);;
131 strncpy(filename,AliDAQ::DdlFileName("PMD",iddl),80);
133 outfile = new AliFstream(filename);
135 // Write the Dummy Data Header into the file
136 Int_t bHPosition = outfile->Tellp();
137 outfile->WriteBuffer((char*)(&header),sizeof(header));
139 for (Int_t ibus = 0; ibus < kbusSize; ibus++)
141 contentsBus[ibus] = 0;
142 for (Int_t ich = 0; ich < kSize; ich++)
144 busPatch[ibus][ich] = 0;
148 for(Int_t ium = 0; ium < 12; ium++)
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
154 mmodule = modulenoddl[ium];
155 if(mmodule == -1) continue;
156 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];
175 // Add the patch Bus header to the DSP contents
176 dsp[i] += 4*dspBus[i];
179 Int_t dspBlockARDL = 0;
180 Int_t dspBlockBRDL = 0;
184 for (Int_t i = 0; i < 5; i++)
187 Int_t iodd = 2*i + 1;
190 dspBlockARDL += dsp[ieven];
191 remainder = dsp[ieven]%2;
199 dspBlockBRDL += dsp[iodd];
200 remainder = dsp[iodd]%2;
211 // Start writing the DDL file
213 AliPMDBlockHeader blHeader;
214 AliPMDDspHeader dspHeader;
215 AliPMDPatchBusHeader pbusHeader;
217 const Int_t kblHLen = blHeader.GetHeaderLength();
218 const Int_t kdspHLen = dspHeader.GetHeaderLength();
219 const Int_t kpbusHLen = pbusHeader.GetHeaderLength();
222 UInt_t dspBlockHeaderWord[8];
223 UInt_t dspHeaderWord[10];
224 UInt_t patchBusHeaderWord[4];
226 UInt_t ddlEndWord[2] = {0xDEADFACE, 0xDEADFACE};
231 for (Int_t iblock = 0; iblock < 2; iblock++)
235 for (Int_t i=0; i<kblHLen; i++)
237 dspBlockHeaderWord[i] = 0;
241 dspBlockHeaderWord[1] = (UInt_t) (dspBlockARDL + kblHLen);
242 dspBlockHeaderWord[2] = (UInt_t) dspBlockARDL;
244 else if (iblock == 1)
246 dspBlockHeaderWord[1] = (UInt_t) (dspBlockBRDL + kblHLen);
247 dspBlockHeaderWord[2] = (UInt_t) dspBlockBRDL;
250 outfile->WriteBuffer((char*)dspBlockHeaderWord,kblHLen*sizeof(UInt_t));
260 else if (iblock == 1)
269 for (Int_t idsp = 0; idsp < 5; idsp++)
276 dspRDL = (UInt_t) dsp[dspno];
278 else if (iblock == 1)
281 dspRDL = (UInt_t) dsp[dspno];
284 for (Int_t i=0; i<kdspHLen; i++)
286 dspHeaderWord[i] = 0;
288 remainder = dspRDL%2;
289 if (remainder == 1) dspRDL++;
291 dspHeaderWord[1] = dspRDL + kdspHLen;
292 dspHeaderWord[2] = dspRDL;
293 dspHeaderWord[3] = dspno;
294 if (remainder == 1) dspHeaderWord[8] = 1; // setting the padding word
297 outfile->WriteBuffer((char*)dspHeaderWord,kdspHLen*sizeof(UInt_t));
299 for (Int_t ibus = 0; ibus < 5; ibus++)
303 Int_t busno = iskip[idsp] + ibus + 1;
304 Int_t patchbusRDL = contentsBus[busno];
308 patchBusHeaderWord[0] = 0;
309 patchBusHeaderWord[1] = (UInt_t) (patchbusRDL + kpbusHLen);
310 patchBusHeaderWord[2] = (UInt_t) patchbusRDL;
311 patchBusHeaderWord[3] = (UInt_t) busno;
313 else if (patchbusRDL == 0)
315 patchBusHeaderWord[0] = 0;
316 patchBusHeaderWord[1] = (UInt_t) kpbusHLen;
317 patchBusHeaderWord[2] = (UInt_t) 0;
318 patchBusHeaderWord[3] = (UInt_t) busno;
322 outfile->WriteBuffer((char*)patchBusHeaderWord,4*sizeof(UInt_t));
324 bknJunk += patchbusRDL;
326 for (Int_t iword = 0; iword < patchbusRDL; iword++)
328 buffer[iword] = busPatch[busno][iword];
331 outfile->WriteBuffer((char*)buffer,patchbusRDL*sizeof(UInt_t));
333 } // End of patch bus loop
336 // Adding a padding word if the total words odd
339 UInt_t paddingWord = dspHeader.GetDefaultPaddingWord();
340 outfile->WriteBuffer((char*)(&paddingWord),sizeof(UInt_t));
345 // Write two extra word at the end of each DDL file
346 outfile->WriteBuffer((char*)ddlEndWord,2*sizeof(UInt_t));
348 // Write real data header
349 // take the pointer to the beginning of the data header
350 // write the total number of words per ddl and bring the
351 // pointer to the current file position and close it
352 UInt_t cFPosition = outfile->Tellp();
353 sizeRawData = cFPosition - bHPosition - sizeof(header);
355 header.fSize = cFPosition - bHPosition;
356 header.SetAttribute(0); // valid data
357 outfile->Seekp(bHPosition);
358 outfile->WriteBuffer((char*)(&header),sizeof(header));
359 outfile->Seekp(cFPosition);
365 //____________________________________________________________________________
367 void AliPMDDDLRawData::GetUMDigitsData(TTree *treeD, Int_t imodule,
368 Int_t ddlno, Int_t *contentsBus,
369 UInt_t busPatch[][1536])
371 // Retrieves digits data UnitModule by UnitModule
373 const Int_t kMaxBus = 51;
376 UInt_t mcmno = 0, chno = 0;
378 Int_t det = 0, smn = 0, irow = 0, icol = 0;
381 // Int_t totPatchBus = 0, bPatchBus = 0, ePatchBus = 0;
382 // Int_t ibus = 0, totmcm = 0, rows = 0, cols = 0, rowe = 0, cole = 0;
383 // Int_t moduleno = 0;
385 Int_t patchBusNo[kMaxBus], mcmperBus[kMaxBus];
386 Int_t startRowBus[kMaxBus], startColBus[kMaxBus];
387 Int_t endRowBus[kMaxBus], endColBus[kMaxBus];
389 Int_t beginPatchBus = -1;
390 Int_t endPatchBus = -1;
391 for(Int_t i = 0; i < kMaxBus; i++)
401 // Fetch the DDL mapping info from the mapping database
403 DdlMapping(ddlno, imodule, beginPatchBus, endPatchBus,
404 patchBusNo, mcmperBus, startRowBus, endRowBus,
405 startColBus, endColBus);
407 // Read if some chains are off from the ddlinfo database
409 Int_t srowoff1[2][24], erowoff1[2][24];
410 Int_t scoloff1[2][24], ecoloff1[2][24];
411 Int_t srowoff2[2][24], erowoff2[2][24];
412 Int_t scoloff2[2][24], ecoloff2[2][24];
414 for (Int_t idet = 0; idet < 2; idet++)
416 for (Int_t im = 0; im < 24; im++)
418 srowoff1[idet][im] = fDdlinfo->GetStartRowA(idet,im);
419 erowoff1[idet][im] = fDdlinfo->GetEndRowA(idet,im);
420 scoloff1[idet][im] = fDdlinfo->GetStartColA(idet,im);
421 ecoloff1[idet][im] = fDdlinfo->GetEndColA(idet,im);
422 srowoff2[idet][im] = fDdlinfo->GetStartRowB(idet,im);
423 erowoff2[idet][im] = fDdlinfo->GetEndRowB(idet,im);
424 scoloff2[idet][im] = fDdlinfo->GetStartColB(idet,im);
425 ecoloff2[idet][im] = fDdlinfo->GetEndColB(idet,im);
429 treeD->GetEntry(imodule);
430 Int_t nentries = fDigits->GetLast();
431 Int_t totword = nentries+1;
433 AliPMDdigit *pmddigit = 0x0;
435 for (Int_t ient = 0; ient < totword; ient++)
437 pmddigit = (AliPMDdigit*)fDigits->UncheckedAt(ient);
439 det = pmddigit->GetDetector();
440 smn = pmddigit->GetSMNumber();
441 irow = pmddigit->GetRow();
442 icol = pmddigit->GetColumn();
443 adc = (UInt_t) pmddigit->GetADC();
445 TransformS2H(smn,irow,icol);
447 // remove the non-existence channels
449 //printf("%d %d %d %d\n",det,smn,irow,icol);
450 //printf("--- %d %d %d %d\n",srowoff[det][smn],erowoff[det][smn],
451 // scoloff[det][smn],ecoloff[det][smn]);
453 if (irow >= srowoff1[det][smn] && irow <= erowoff1[det][smn])
455 if (icol >= scoloff1[det][smn] && icol <= ecoloff1[det][smn])
460 if (irow >= srowoff2[det][smn] && irow <= erowoff2[det][smn])
462 if (icol >= scoloff2[det][smn] && icol <= ecoloff2[det][smn])
469 GetMCMCh(imodule, irow, icol, beginPatchBus, endPatchBus,
470 mcmperBus, startRowBus, startColBus,
471 endRowBus, endColBus, busno, mcmno, chno);
474 AliBitPacking::PackWord(adc,baseword,0,11);
475 AliBitPacking::PackWord(chno,baseword,12,17);
476 AliBitPacking::PackWord(mcmno,baseword,18,28);
477 AliBitPacking::PackWord(0,baseword,29,30);
478 parity = ComputeParity(baseword); // generate the parity bit
479 AliBitPacking::PackWord(parity,baseword,31,31);
481 Int_t jj = contentsBus[busno];
482 busPatch[busno][jj] = baseword;
484 contentsBus[busno]++;
489 //____________________________________________________________________________
490 void AliPMDDDLRawData::TransformS2H(Int_t smn, Int_t &irow, Int_t &icol)
492 // Does the Software to Hardware coordinate transformation
498 // First in digits we have all dimension 48x96
499 // Transform into the realistic one, i.e, For SM 0&1 96(row)x48(col)
500 // and for SM 2&3 48(row)x96(col)
507 else if( smn >= 12 && smn < 24)
519 //____________________________________________________________________________
521 void AliPMDDDLRawData::GetMCMCh(Int_t imodule, Int_t row, Int_t col,
522 Int_t beginPatchBus, Int_t endPatchBus,
524 Int_t *startRowBus, Int_t *startColBus,
525 Int_t *endRowBus, Int_t *endColBus,
526 Int_t & busno, UInt_t &mcmno, UInt_t &chno)
528 // This converts row col to hardware channel number
529 // This is the final version of mapping supplied by Mriganka
533 static const UInt_t kChDdl01[16][4] = { {6, 4, 5, 7},
551 static const UInt_t kChDdl23[16][4] = { {57, 59, 58, 56},
569 static const UInt_t kChDdl41[16][4] = { {56, 58, 59, 57},
587 static const UInt_t kChDdl42[16][4] = { {7, 5, 4, 6},
605 static const UInt_t kChDdl51[16][4] = { {7, 5, 4, 6},
624 static const UInt_t kChDdl52[16][4] = { {56, 58, 59, 57},
642 for (Int_t i = 0; i < 16; i++)
644 for (Int_t j = 0; j < 4; j++)
647 if(imodule < 6) iCh[i][j] = kChDdl01[i][j];
648 if(imodule >= 6 && imodule <= 11) iCh[i][j] = kChDdl01[i][j];
649 if(imodule >= 12 && imodule <= 17) iCh[i][j] = kChDdl23[i][j];
650 if(imodule >= 18 && imodule <= 23) iCh[i][j] = kChDdl23[i][j];
652 if(imodule >= 24 && imodule <= 29) iCh[i][j] = kChDdl41[i][j];
653 if(imodule >= 42 && imodule <= 47) iCh[i][j] = kChDdl42[i][j];
654 if(imodule >= 36 && imodule <= 41) iCh[i][j] = kChDdl51[i][j];
655 if(imodule >= 30 && imodule <= 35) iCh[i][j] = kChDdl52[i][j];
661 Int_t irownew = row%16;
662 Int_t icolnew = col%4;
664 chno = iCh[irownew][icolnew];
667 for (Int_t ibus = beginPatchBus; ibus <= endPatchBus; ibus++)
669 Int_t srow = startRowBus[ibus];
670 Int_t erow = endRowBus[ibus];
671 Int_t scol = startColBus[ibus];
672 Int_t ecol = endColBus[ibus];
673 Int_t tmcm = mcmperBus[ibus];
675 if ((row >= srow && row <= erow) && (col >= scol && col <= ecol))
679 // Find out the MCM Number
682 if (imodule < 6) mcmno = (col-scol)/4 + 1;
683 if (imodule >= 6 && imodule < 12) mcmno = (col-scol)/4 + 1;
685 if (imodule >= 12 && imodule < 18)
687 icolnew = (col - scol)/4;
688 mcmno = tmcm - icolnew;
690 if (imodule >= 18 && imodule < 24)
692 icolnew = (col - scol)/4;
693 mcmno = tmcm - icolnew;
697 if (imodule >= 24 && imodule < 30)
701 Int_t rowdiff = endRowBus[ibus] - startRowBus[ibus];
704 Int_t midrow = srow + 16;
705 if(row >= srow && row < midrow)
707 mcmno = 12 + (col-scol)/4 + 1;
709 else if(row >= midrow && row <= erow)
712 mcmno = (col-scol)/4 + 1;
715 else if (rowdiff < 16)
717 mcmno = (col-scol)/4 + 1;
721 if (imodule >= 42 && imodule < 48)
723 Int_t rowdiff = endRowBus[ibus] - startRowBus[ibus];
726 Int_t midrow = srow + 16;
727 if (row >= midrow && row <= erow)
729 mcmno = 12 + (ecol -col)/4 + 1;
731 else if (row >= srow && row < midrow)
733 mcmno = (ecol - col)/4 + 1;
736 else if (rowdiff < 16)
738 mcmno = (ecol - col)/4 + 1;
743 if (imodule >= 30 && imodule < 36)
745 // CPV plane, SU Mod = 1, 2 : ddl = 5
748 Int_t rowdiff = endRowBus[ibus] - startRowBus[ibus];
751 Int_t midrow = srow + 16;
752 if(row >= srow && row < midrow)
754 mcmno = 12 + (col-scol)/4 + 1;
756 else if(row >= midrow && row <= erow)
758 mcmno = (col-scol)/4 + 1;
761 else if(rowdiff < 16)
763 mcmno = (col-scol)/4 + 1;
767 if (imodule >= 36 && imodule < 42)
769 Int_t rowdiff = endRowBus[ibus] - startRowBus[ibus];
772 Int_t midrow = srow + 16;
773 if (row >= midrow && row <= erow)
775 mcmno = 12 + (ecol - col)/4 + 1;
777 else if (row >= srow && row < midrow)
779 mcmno = (ecol - col)/4 + 1;
782 else if (rowdiff < 16)
784 mcmno = (ecol - col)/4 + 1;
792 //____________________________________________________________________________
794 Int_t AliPMDDDLRawData::ComputeParity(UInt_t baseword)
796 // Generate the parity bit
799 for(Int_t j=0; j<29; j++)
801 if (baseword & 0x01 ) count++;
804 Int_t parity = count%2;
807 //____________________________________________________________________________
808 void AliPMDDDLRawData::DdlMapping(Int_t iddl, Int_t imodule,
809 Int_t &beginPatchBus, Int_t &endPatchBus,
810 Int_t patchBusNo[], Int_t mcmperBus[],
811 Int_t startRowBus[], Int_t endRowBus[],
812 Int_t startColBus[], Int_t endColBus[])
814 // DDL Mapping fetching from mapping database
816 beginPatchBus = fMapData->GetBeginPatchBus(iddl,imodule);
817 endPatchBus = fMapData->GetEndPatchBus(iddl,imodule);
819 for(Int_t ibus = beginPatchBus; ibus < endPatchBus+1; ibus++)
821 patchBusNo[ibus] = ibus;
822 mcmperBus[ibus] = fMapData->GetMcmperBus(iddl,ibus);
823 startRowBus[ibus] = fMapData->GetStartRowBus(iddl,ibus);
824 startColBus[ibus] = fMapData->GetStartColBus(iddl,ibus);
825 endRowBus[ibus] = fMapData->GetEndRowBus(iddl,ibus);
826 endColBus[ibus] = fMapData->GetEndColBus(iddl,ibus);
830 //____________________________________________________________________________
832 AliPMDddlinfoData* AliPMDDDLRawData::GetDdlinfoData() const
834 // The run number will be centralized in AliCDBManager,
835 // you don't need to set it here!
836 AliCDBEntry *entry = AliCDBManager::Instance()->Get("PMD/Calib/Ddlinfo");
838 if(!entry) AliFatal("ddlinfo object retrieval failed!");
840 AliPMDddlinfoData *ddlinfo = 0;
841 if (entry) ddlinfo = (AliPMDddlinfoData*) entry->GetObject();
843 if (!ddlinfo) AliFatal("No ddl info data from database !");
847 //____________________________________________________________________________
848 AliPMDMappingData* AliPMDDDLRawData::GetMappingData() const
850 // The run number will be centralized in AliCDBManager,
851 // you don't need to set it here!
852 AliCDBEntry *entry = AliCDBManager::Instance()->Get("PMD/Calib/Mapping");
854 if(!entry) AliFatal("Mapping object retrieval failed!");
856 AliPMDMappingData *mapda = 0;
857 if (entry) mapda = (AliPMDMappingData*) entry->GetObject();
859 if (!mapda) AliFatal("No mapping data from database !");
863 //____________________________________________________________________________