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>
23 #include "AliRawDataHeader.h"
24 #include "AliBitPacking.h"
25 #include "AliPMDdigit.h"
26 #include "AliPMDRawStream.h"
27 #include "AliPMDDDLRawData.h"
30 ClassImp(AliPMDDDLRawData)
32 AliPMDDDLRawData::AliPMDDDLRawData():
33 fDigits(new TClonesArray("AliPMDdigit", 1000))
35 // Default Constructor
39 //____________________________________________________________________________
41 AliPMDDDLRawData::~AliPMDDDLRawData()
48 //____________________________________________________________________________
49 void AliPMDDDLRawData::WritePMDRawData(TTree *treeD)
51 // write digits into raw data format
55 TBranch *branch = treeD->GetBranch("PMDDigit");
58 AliError("PMD Digit branch not found");
61 branch->SetAddress(&fDigits);
63 Int_t nmodules = (Int_t) treeD->GetEntries();
64 AliDebug(1,Form("Number of modules inside treeD = %d",nmodules));
66 const Int_t kSize = 4608;
68 Int_t modulePerDDL = 0;
71 AliRawDataHeader header;
72 UInt_t sizeRawData = 0;
80 for(Int_t iddl = 0; iddl < kDDL; iddl++)
82 sprintf(filename,"PMD_%d.ddl",AliPMDRawStream::kDDLOffset+iddl);
84 outfile.open(filename,ios::binary);
86 outfile.open(filename);
106 // Write the Dummy Data Header into the file
107 Int_t bHPosition = outfile.tellp();
108 outfile.write((char*)(&header),sizeof(header));
111 for(Int_t ium = 0; ium < modulePerDDL; ium++)
114 if (iddl == 4 && ium == 6) mmodule = 36;
115 if (iddl == 5 && ium == 6) mmodule = 42;
117 for (Int_t i = 0; i < kSize; i++)
121 // Extract energy deposition per cell and pack it
122 // in a 32-bit word and returns all the total words
123 // per one unit-module
125 GetUMDigitsData(treeD, mmodule, ium, iddl, totword, buffer);
127 outfile.write((char*)buffer,totword*sizeof(UInt_t));
133 // Write real data header
134 // take the pointer to the beginning of the data header
135 // write the total number of words per ddl and bring the
136 // pointer to the current file position and close it
137 UInt_t cFPosition = outfile.tellp();
138 sizeRawData = cFPosition - bHPosition - sizeof(header);
139 header.fSize = cFPosition - bHPosition;
140 header.SetAttribute(0); // valid data
141 outfile.seekp(bHPosition);
142 outfile.write((char*)(&header),sizeof(header));
143 outfile.seekp(cFPosition);
150 //____________________________________________________________________________
151 void AliPMDDDLRawData::GetUMDigitsData(TTree *treeD, Int_t imodule, Int_t ium,
152 Int_t ddlno, Int_t & totword,
155 // Retrives digits data UnitModule by UnitModule
159 Int_t det, smn, irow, icol;
161 treeD->GetEntry(imodule);
162 Int_t nentries = fDigits->GetLast();
163 totword = nentries+1;
165 for (Int_t ient = 0; ient < totword; ient++)
167 fPMDdigit = (AliPMDdigit*)fDigits->UncheckedAt(ient);
169 det = fPMDdigit->GetDetector();
170 smn = fPMDdigit->GetSMNumber();
171 irow = fPMDdigit->GetRow();
172 icol = fPMDdigit->GetColumn();
173 adc = (UInt_t) fPMDdigit->GetADC();
175 TransformS2H(smn,irow,icol);
176 GetMCMCh(ddlno, ium, irow, icol, mcmno, chno);
179 AliBitPacking::PackWord(adc,baseword,0,11);
180 AliBitPacking::PackWord(chno,baseword,12,17);
181 AliBitPacking::PackWord(mcmno,baseword,18,28);
182 AliBitPacking::PackWord(0,baseword,29,30);
183 AliBitPacking::PackWord(1,baseword,31,31);
185 buffer[ient] = baseword;
190 //____________________________________________________________________________
191 void AliPMDDDLRawData::TransformS2H(Int_t smn, Int_t &irow, Int_t &icol)
193 // Does the Software to Hardware coordinate transformation
200 // First in digits we have all dimension 48x96
201 // Transform into the realistic one, i.e, For SM 1&2 96x48
202 // and for SM 3&4 48x96
209 else if( smn >= 12 && smn < 24)
214 // This is the transformation of Geant (0,0) to the Hardware (0,0)
215 // which is always at the top left corner. This may change in future.
216 // Then accordingly we have to transform it.
219 irownew = 95 - irownew1;
222 else if(smn >= 6 && smn < 12)
225 icolnew = 47 - icolnew1;
227 else if(smn >= 12 && smn < 18)
229 irownew = 47 - irownew1;
232 else if(smn >= 18 && smn < 24)
235 icolnew = 95 - icolnew1;
244 //____________________________________________________________________________
246 void AliPMDDDLRawData::GetMCMCh(Int_t ddlno, Int_t um,
247 Int_t row, Int_t col,
248 UInt_t &mcmno, UInt_t &chno)
250 // This part will be modified once the final track layout on the PCB is
251 // designed. This will only change the coordinate of the individual cell
253 static const UInt_t kCh[16][4] = { {12, 13, 18, 19},
270 if (ddlno == 0 || ddlno == 1)
272 // PRE plane, SU Mod = 1, 2
273 Int_t irownew = row%16;
274 Int_t icolnew = col%4;
275 Int_t irowdiv = row/16;
276 Int_t icoldiv = col/4;
278 mcmno = 72*um + 12*irowdiv + icoldiv;
279 chno = kCh[irownew][icolnew];
281 else if (ddlno == 2 || ddlno == 3)
283 // PRE plane, SU Mod = 3, 4
284 Int_t irownew = row%16;
285 Int_t icolnew = col%4;
286 Int_t irowdiv = row/16;
287 Int_t icoldiv = col/4;
289 mcmno = 72*um + 24*irowdiv + icoldiv;
290 chno = kCh[irownew][icolnew];
292 else if (ddlno == 4 || ddlno == 5)
294 // CPV plane, SU Mod = 1, 3 : ddl = 4
295 // CPV plane, SU Mod = 2, 4 : ddl = 5
297 Int_t irownew = row%16;
298 Int_t icolnew = col%4;
299 Int_t irowdiv = row/16;
300 Int_t icoldiv = col/4;
304 mcmno = 72*um + 12*irowdiv + icoldiv;
308 mcmno = 72*um + 24*irowdiv + icoldiv;
310 chno = kCh[irownew][icolnew];
314 //____________________________________________________________________________