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>
22 #include "AliPMDdigit.h"
23 #include "AliPMDDDLRawData.h"
26 ClassImp(AliPMDDDLRawData)
28 AliPMDDDLRawData::AliPMDDDLRawData():
29 fDigits(new TClonesArray("AliPMDdigit", 1000))
31 // Default Constructor
34 //____________________________________________________________________________
36 AliPMDDDLRawData::~AliPMDDDLRawData()
43 //____________________________________________________________________________
44 void AliPMDDDLRawData::WritePMDRawData(TTree *treeD, Int_t evtno)
46 // write digits into raw data format
50 TBranch *branch = treeD->GetBranch("PMDDigit");
52 branch->SetAddress(&fDigits);
54 Int_t nmodules = (Int_t) treeD->GetEntries();
55 // cout << " nmodules = " << nmodules << endl;
57 const Int_t kSize = 4608;
59 Int_t modulePerDDL = 0;
62 UInt_t sizeRawData = 0;
63 UInt_t magicWord = 0x123456;
64 UInt_t detectorID = 0;
68 UInt_t mHSize = 3*sizeof(UInt_t);
77 for(Int_t iddl = 0; iddl < kDDL; iddl++)
79 sprintf(filename,"Ev%dPMDddl%d.dat",evtno,iddl);
81 outfile.open(filename,ios::binary);
83 outfile.open(filename);
105 miniHeader[0] = sizeRawData;
106 PackWord(0,23,magicWord,miniHeader[1]);
107 PackWord(24,31,detectorID,miniHeader[1]);
109 PackWord(0,15,ddlID,miniHeader[2]);
110 PackWord(16,23,flag,miniHeader[2]);
111 PackWord(24,31,version,miniHeader[2]);
114 // Write the Dummy Mini Header into the file
115 Int_t bHPosition = outfile.tellp();
116 outfile.write((char*)(miniHeader),mHSize);
119 for(Int_t ium = 0; ium < modulePerDDL; ium++)
122 for (Int_t i = 0; i < kSize; i++)
126 // Extract energy deposition per cell and pack it
127 // in a 32-bit word and returns all the total words
128 // per one unit-module
130 GetUMDigitsData(treeD, mmodule, ium, iddl, totword, buffer);
132 outfile.write((char*)buffer,totword*sizeof(UInt_t));
138 // Write real mini header
139 // take the pointer to the beginning of the mini header
140 // write the total number of words per ddl and bring the
141 // pointer to the current file position and close it
142 UInt_t cFPosition = outfile.tellp();
143 sizeRawData = cFPosition - bHPosition - mHSize;
144 outfile.seekp(bHPosition);
145 outfile.write((char*)(&sizeRawData),sizeof(UInt_t));
146 outfile.seekp(cFPosition);
153 //____________________________________________________________________________
155 void AliPMDDDLRawData::ReadPMDRawData(Int_t evtno)
157 // reads the raw data
162 const Int_t kDDL = 6;
165 for(Int_t iddl = 0; iddl < kDDL; iddl++)
167 sprintf(filename,"Ev%dPMDddl%d.dat",evtno,iddl);
169 infile.open(filename,ios::binary);
171 infile.open(filename);
178 UInt_t miniHeader[3];
180 UInt_t mHSize = 3*sizeof(UInt_t);
182 infile.read((char*)(miniHeader), mHSize);
184 sizeRawData = miniHeader[0];
185 Int_t totword = sizeRawData/4;
191 for (Int_t iword = 0; iword < totword; iword++)
193 infile.read((char*)(&baseword),sizeof(UInt_t));
195 UnpackWord(0,11,adc,baseword);
196 UnpackWord(12,17,chno,baseword);
197 UnpackWord(18,28,mcmno,baseword);
199 GetRowCol(iddl, mcmno, chno, ium, irow, icol);
202 imodule = iddl*6 + ium;
213 fprintf(fpw3,"%d %d %d %d %d %d\n",imodule,irow,icol,mcmno,chno,adc);
219 //____________________________________________________________________________
220 void AliPMDDDLRawData::GetUMDigitsData(TTree *treeD, Int_t imodule, Int_t ium,
221 Int_t ddlno, Int_t & totword,
225 UInt_t dataword, baseword;
231 Int_t det, smn, irow, icol;
234 treeD->GetEntry(imodule);
235 Int_t nentries = fDigits->GetLast();
236 totword = nentries+1;
238 for (Int_t ient = 0; ient < totword; ient++)
240 fPMDdigit = (AliPMDdigit*)fDigits->UncheckedAt(ient);
242 det = fPMDdigit->GetDetector();
243 smn = fPMDdigit->GetSMNumber();
244 irow = fPMDdigit->GetRow();
245 icol = fPMDdigit->GetColumn();
246 adc = (UInt_t) fPMDdigit->GetADC();
253 else if( smn >= 12 && smn < 24)
261 GetMCMCh(ddlno, ium, irownew, icolnew, mcmno, chno);
263 fprintf(fpw2,"%d %d %d %d %d %d\n",imodule,irownew,icolnew,mcmno,chno,adc);
267 PackWord(0, 11,dataword,baseword);
269 PackWord(12,17,dataword,baseword);
271 PackWord(18,28,dataword,baseword);
273 PackWord(29,29,dataword,baseword);
275 PackWord(30,30,dataword,baseword);
277 PackWord(31,31,dataword,baseword);
279 buffer[ient] = baseword;
283 // fprintf(fpw1,"%d %d %d %d %d\n",det,smn,irow,icol,adc);
289 //____________________________________________________________________________
290 void AliPMDDDLRawData::GetMCMCh(Int_t ddlno, Int_t um,
291 Int_t row, Int_t col,
292 UInt_t &mcmno, UInt_t &chno)
294 // This part will be modified once the final track layout on the PCB is
295 // designed. This will only change the coordinate of the individual cell
296 UInt_t ch[16][4] = { {3, 2, 29, 28},
313 if (ddlno == 0 || ddlno == 1)
315 // PRE plane, SU Mod = 1, 2
316 Int_t irownew = row%16;
317 Int_t icolnew = col%4;
318 Int_t irowdiv = row/16;
319 Int_t icoldiv = col/4;
321 mcmno = 72*um + 6*icoldiv + irowdiv;
322 chno = ch[irownew][icolnew];
324 else if (ddlno == 2 || ddlno == 3)
326 // PRE plane, SU Mod = 3, 4
327 Int_t irownew = row%16;
328 Int_t icolnew = col%4;
329 Int_t irowdiv = row/16;
330 Int_t icoldiv = col/4;
332 mcmno = 72*um + 3*icoldiv + irowdiv;
333 chno = ch[irownew][icolnew];
337 // CPV plane, SU Mod = 1, 2
338 Int_t irownew = row%16;
339 Int_t icolnew = col%4;
340 Int_t irowdiv = row/16;
341 Int_t icoldiv = col/4;
343 mcmno = 72*um + 6*icoldiv + irowdiv;
344 chno = ch[irownew][icolnew];
348 // CPV plane, SU Mod = 3, 4
349 Int_t irownew = row%16;
350 Int_t icolnew = col%4;
351 Int_t irowdiv = row/16;
352 Int_t icoldiv = col/4;
354 mcmno = 72*um + 3*icoldiv + irowdiv;
355 chno = ch[irownew][icolnew];
359 //____________________________________________________________________________
361 void AliPMDDDLRawData::GetRowCol(Int_t ddlno, UInt_t mcmno, UInt_t chno,
362 Int_t &um, Int_t &row, Int_t &col)
364 UInt_t ch[64] = { 9, 5, 1, 0, 13, 4, 8, 12,
365 16, 17, 20, 24, 28, 29, 21, 25,
366 22, 26, 30, 31, 18, 27, 23, 19,
367 15, 14, 11, 7, 3, 2, 10, 6,
368 41, 37, 33, 32, 45, 36, 40, 44,
369 48, 49, 52, 56, 60, 61, 53, 57,
370 54, 58, 62, 63, 50, 59, 55, 51,
371 47, 46, 43, 39, 35, 34, 42, 38 };
373 if (ddlno == 1 || ddlno == 2)
376 Int_t mcmnonew = mcmno - 72*um;
377 Int_t rowcol = ch[chno];
378 Int_t irownew = rowcol/4;
379 Int_t icolnew = rowcol%4;
381 Int_t remmcm = mcmnonew%6;
382 Int_t divmcm = mcmnonew/6;
384 row = 16*remmcm + irownew;
385 col = 4*divmcm + icolnew;
387 else if (ddlno == 3 || ddlno == 4)
390 Int_t mcmnonew = mcmno - 72*um;
391 Int_t rowcol = ch[chno];
392 Int_t irownew = rowcol/4;
393 Int_t icolnew = rowcol%4;
395 Int_t remmcm = mcmnonew%3;
396 Int_t divmcm = mcmnonew/3;
398 row = 16*remmcm + irownew;
399 col = 4*divmcm + icolnew;
404 Int_t mcmnonew = mcmno - 72*um;
405 Int_t rowcol = ch[chno];
406 Int_t irownew = rowcol/4;
407 Int_t icolnew = rowcol%4;
409 Int_t remmcm = mcmnonew%6;
410 Int_t divmcm = mcmnonew/6;
412 row = 16*remmcm + irownew;
413 col = 4*divmcm + icolnew;
418 Int_t mcmnonew = mcmno - 72*um;
419 Int_t rowcol = ch[chno];
420 Int_t irownew = rowcol/4;
421 Int_t icolnew = rowcol%4;
423 Int_t remmcm = mcmnonew%3;
424 Int_t divmcm = mcmnonew/3;
426 row = 16*remmcm + irownew;
427 col = 4*divmcm + icolnew;
432 //____________________________________________________________________________
434 void AliPMDDDLRawData::PackWord(UInt_t startbit, UInt_t stopbit,
435 UInt_t dataword, UInt_t &packedword)
437 UInt_t bitLength = stopbit - startbit + 1;
438 UInt_t bitContent = (UInt_t) (TMath::Power(2,bitLength) - 1);
439 if(bitContent < dataword)
441 cout << " *** ERROR *** bitContent is less than the dataword" << endl;
444 UInt_t packedBits = 0;
446 packedBits = (UInt_t) (TMath::Log(packedword)/TMath::Log(2));
449 if (packedBits <= stopbit)
451 counter = 31 - stopbit;
455 counter = 31 - packedBits;
457 UInt_t dummyword = 0xFFFFFFFF;
458 dummyword >>= counter;
459 UInt_t lword = dataword << startbit;
460 UInt_t nword = lword | packedword;
461 packedword = dummyword & nword;
465 //____________________________________________________________________________
466 void AliPMDDDLRawData::UnpackWord(UInt_t startbit, UInt_t stopbit,
467 UInt_t &dataword, UInt_t packedword)
469 UInt_t bitLength = stopbit - startbit + 1;
470 UInt_t bitContent = (UInt_t) (TMath::Power(2,bitLength) - 1);
471 bitContent <<= startbit;
472 dataword = packedword & bitContent;
473 dataword >>= startbit;
476 //____________________________________________________________________________