CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DTTFFEDSim.cc
Go to the documentation of this file.
1 //-------------------------------------------------
2 //
3 // Class: DTTFFEDSim
4 //
5 // L1 DT Track Finder Digi-to-Raw
6 //
7 //
8 //
9 // Author :
10 // J. Troconiz UAM Madrid
11 //
12 //--------------------------------------------------
13 
15 
18 
20 
21 #include <iostream>
22 
23 using namespace std;
24 
25 DTTFFEDSim::DTTFFEDSim(const edm::ParameterSet& pset) : eventNum(0) {
26 
27  produces<FEDRawDataCollection>();
28 
29  DTDigiInputTag = pset.getParameter<edm::InputTag>("DTDigi_Source");
30  DTPHTFInputTag = pset.getParameter<edm::InputTag>("DTTracks_Source");
31 
32  ChPh_tok = consumes<L1MuDTChambPhContainer>(DTDigiInputTag);
33  ChTh_tok = consumes<L1MuDTChambThContainer>(DTDigiInputTag);
34  Trk_tok = consumes<L1MuDTTrackContainer>(DTPHTFInputTag);
35 
36 }
37 
39 
41 
43 
44  if (!fillRawData(e, data)) return;
45 
46  auto_ptr<FEDRawDataCollection> fed_product(new FEDRawDataCollection(data));
47 
48  e.put(fed_product);
49 
50 }
51 
54 
55  eventNum = e.id().event();
56 
57  int lines = 2;
58 
60  e.getByToken(ChPh_tok,phtrig);
61  lines += phtrig->bxSize(-1, 1);
62 
64  e.getByToken(ChTh_tok,thtrig);
65  lines += thtrig->bxSize(-1, 1);
66 
68  e.getByToken(Trk_tok,trtrig);
69  lines += trtrig->bxSize(-1, 1)*3;
70 
71  FEDRawData& dttfdata = data.FEDData(0x30C);
72  dttfdata.resize(lines*8); // size in bytes
73  unsigned char* LineFED=dttfdata.data();
74 
75  int* dataWord1 = new int;
76  int* dataWord2 = new int;
77 
78  //--> Header
79 
80  *dataWord1 = 0x50000000
81  + (eventNum&0xFFFFFF);
82  *dataWord2 = 0x00030C00;
83 
84  int newCRC = 0xFFFF;
85  calcCRC(*dataWord1, *dataWord2, newCRC);
86  *((int*)LineFED)=*dataWord2;
87  LineFED+=4;
88  *((int*)LineFED)=*dataWord1;
89 
90  //--> DTTF data
91 
92  int TS1Id[4], TS2Id[4]; // word identifier for TS #1,#2 for stations
93  TS1Id[0] = 0x0E;
94  TS2Id[0] = 0x1E;
95  TS1Id[1] = 0x2E;
96  TS2Id[1] = 0x3E;
97  TS1Id[3] = 0x4E;
98  TS2Id[3] = 0x5E;
99  TS1Id[2] = 0x8FFF8;
100  TS2Id[2] = 0x9FFF8;
101 
102  //Input
104 
105  for ( tsphi = phtrig->getContainer()->begin();
106  tsphi != phtrig->getContainer()->end();
107  tsphi++ ) {
108  if ( tsphi->code() != 7 ) {
109 
110  int wheelID = tsphi->whNum()+1;
111  if ( wheelID <= 0 ) wheelID -= 2;
112  int stationID = tsphi->stNum()-1;
113  int is2nd = tsphi->Ts2Tag();
114 
115  int channelNr = channel(wheelID, tsphi->scNum(), tsphi->bxNum()-is2nd);
116  if ( channelNr == 255 ) continue;
117  int TSId = ( is2nd == 0 ) ? TS1Id[stationID] : TS2Id[stationID];
118 
119  *dataWord1 = ((channelNr&0xFF)<<24)
120  + 0x00FFFFFF;
121 
122  if ( stationID != 2 ){
123  *dataWord2 = ( (TSId&0x0FF)<<24)
124  + (~(tsphi->code()+1)&0x007)
125  + ( (~tsphi->phiB()&0x3FF)<<3)
126  + ( (~tsphi->phi()&0xFFF)<<13);
127  }
128  else {
129  *dataWord2 = ( (TSId&0xFFFFF)<<12)
130  + (~(tsphi->code()+1)&0x00007)
131  + ( (~tsphi->phi()&0x00FFF)<<3);
132  }
133 
134  calcCRC(*dataWord1, *dataWord2, newCRC);
135  LineFED+=4;
136  *((int*)LineFED)=*dataWord2;
137  LineFED+=4;
138  *((int*)LineFED)=*dataWord1;
139  }
140  }
141  //Input
142 
143  //Input
145 
146  for ( tsthe = thtrig->getContainer()->begin();
147  tsthe != thtrig->getContainer()->end();
148  tsthe++ ) {
149 
150  int wheelTh = tsthe->whNum();
151  int sectorID = tsthe->scNum();
152 
153  int channelNr = channel(0, sectorID, tsthe->bxNum());
154  if ( channelNr == 255 ) continue;
155  int TSId = wheelTh+2;
156 
157  *dataWord1 = ((channelNr&0xFF)<<24)
158  + 0x00FFFFFF;
159 
160  *dataWord2 = ((TSId&0x07)<<28)
161  + 0x0FFFFFFF;
162 
163  int stationID = tsthe->stNum()-1;
164  for ( int bti = 0; bti < 7; bti++ )
165  if ( wheelTh == -2 || wheelTh == -1 ||
166  ( wheelTh == 0 && (sectorID == 0 || sectorID == 3 || sectorID == 4 || sectorID == 7 || sectorID == 8 || sectorID == 11) ) )
167  *dataWord2 -= (tsthe->position(bti)&0x1)<<(stationID*7+bti);
168  else
169  *dataWord2 -= (tsthe->position(6-bti)&0x1)<<(stationID*7+bti);
170 
171  calcCRC(*dataWord1, *dataWord2, newCRC);
172  LineFED+=4;
173  *((int*)LineFED)=*dataWord2;
174  LineFED+=4;
175  *((int*)LineFED)=*dataWord1;
176  }
177  //Input
178 
179  //Output
181 
182  for ( tstrk = trtrig->getContainer()->begin();
183  tstrk != trtrig->getContainer()->end();
184  tstrk++ ) {
185 
186  int channelNr = channel(tstrk->whNum(), tstrk->scNum(), tstrk->bx());
187  if ( channelNr == 255 ) continue;
188  int TSId = ( tstrk->TrkTag() == 0 ) ? 0xAFFF : 0xBFFF;
189 
190  *dataWord1 = ((channelNr&0xFF)<<24)
191  + 0x00FFFFFF;
192 
193  *dataWord2 = ( (TSId&0xFFFF)<<16)
194  + ( tstrk->stNum(4)&0x0000F)
195  + ((tstrk->stNum(3)&0x0000F)<<4)
196  + ((tstrk->stNum(2)&0x0000F)<<8)
197  + ((tstrk->stNum(1)&0x00003)<<12);
198 
199  calcCRC(*dataWord1, *dataWord2, newCRC);
200  LineFED+=4;
201  *((int*)LineFED)=*dataWord2;
202  LineFED+=4;
203  *((int*)LineFED)=*dataWord1;
204 
205  TSId = ( tstrk->TrkTag() == 0 ) ? 0xCFFE : 0xDFFE;
206 
207  *dataWord1 = ((channelNr&0xFF)<<24)
208  + 0x00FFFFFF;
209 
210  *dataWord2 = ( (TSId&0xFFFE)<<16)
211  + ( ~tstrk->quality_packed()&0x0007)
212  + ( (tstrk->phi_packed()&0x00FF)<<3)
213  + ( (~tstrk->charge_packed()&0x0001)<<11)
214  + ( (~tstrk->pt_packed()&0x001F)<<12);
215 
216  calcCRC(*dataWord1, *dataWord2, newCRC);
217  LineFED+=4;
218  *((int*)LineFED)=*dataWord2;
219  LineFED+=4;
220  *((int*)LineFED)=*dataWord1;
221 
222  channelNr = channel(0, tstrk->scNum(), tstrk->bx());
223  if ( channelNr == 255 ) continue;
224  TSId = (tstrk->whNum()+3)<<16;
225  TSId += ( tstrk->whNum() < 0 ) ? 0x8FFFC : 0x7FFFC;
226 
227  *dataWord1 = ((channelNr&0xFF)<<24)
228  + 0x00FFFFFF;
229 
230  *dataWord2 = (TSId&0xFFFFC)<<12;
231 
232  if ( tstrk->TrkTag() == 0 ) {
233  *dataWord2 += 0x3F80
234  + ( tstrk->eta_packed()&0x003F)
235  + ((~tstrk->finehalo_packed()&0x0001)<<6);
236  }
237  else {
238  *dataWord2 += 0x007F
239  + ( ( tstrk->eta_packed()&0x003F)<<7)
240  + ((~tstrk->finehalo_packed()&0x0001)<<13);
241  }
242 
243  calcCRC(*dataWord1, *dataWord2, newCRC);
244  LineFED+=4;
245  *((int*)LineFED)=*dataWord2;
246  LineFED+=4;
247  *((int*)LineFED)=*dataWord1;
248  }
249  //Output
250 
251  //--> Trailer
252 
253  *dataWord1 = 0xA0000000
254  + (lines&0xFFFFFF);
255  *dataWord2 = 0;
256 
257  calcCRC(*dataWord1, *dataWord2&0xFFFF, newCRC);
258 
259  *dataWord2 += (newCRC&0xFFFF)<<16;
260 
261  LineFED+=4;
262  *((int*)LineFED)=*dataWord2;
263  LineFED+=4;
264  *((int*)LineFED)=*dataWord1;
265 
266  delete dataWord1;
267  delete dataWord2;
268  return true;
269 }
270 
271 int DTTFFEDSim::channel( int wheel, int sector, int bx ){
272 
273  // wheel : -3 -2 -1 +1 +2 +3 <=> PHTF's : N2, N1, N0, P0, P1, P2
274  // 0 <=> ETTF
275  // sector : 0 -> 11
276  // bx : -1 -> +1
277 
278  int myChannel = 255;
279 
280  if ( abs(bx) > 1) { return myChannel; }
281  if ( sector < 0 || sector > 11) { return myChannel; }
282  if ( abs(wheel) > 3) { return myChannel; }
283 
284  myChannel = sector*21 + wheel*3 - bx + 10 ;
285 
286  if (myChannel > 125) myChannel += 2;
287 
288  return myChannel;
289 }
290 
291 int DTTFFEDSim::bxNr( int channel ){
292 
293  int myChannel = channel;
294 
295  if (myChannel > 127) myChannel -= 2;
296 
297  if (myChannel < 0 || myChannel > 251 ){ return -999; }
298 
299  int myBx = 1-(myChannel%3);
300 
301  return myBx;
302 }
303 
304 int DTTFFEDSim::sector( int channel ){
305 
306  int myChannel = channel;
307 
308  if (myChannel > 127) myChannel -= 2;
309 
310  if (myChannel < 0 || myChannel > 251 ){ return -999; }
311 
312  return myChannel/21;
313 }
314 
315 int DTTFFEDSim::wheel( int channel ){
316 
317  int myChannel = channel;
318 
319  if (myChannel > 127) myChannel -= 2;
320 
321  if (myChannel < 0 || myChannel > 251 ){ return -999; }
322 
323  int myWheel = ((myChannel%21)/3)-3;
324 
325  return myWheel;
326 }
327 
328 void DTTFFEDSim::calcCRC(int myD1, int myD2, int &myC){
329 
330  int myCRC[16],D[64],C[16];
331 
332  for( int i=0; i < 32; i++ ){ D[i]=(myD2>>i)&0x1; }
333  for( int i=0; i < 32; i++ ){ D[i+32]=(myD1>>i)&0x1; }
334  for( int i=0; i < 16; i++ ){ C[i]=(myC>>i)&0x1; }
335 
336  myCRC[0] = ( D[63] + D[62] + D[61] + D[60] + D[55] + D[54] +
337  D[53] + D[52] + D[51] + D[50] + D[49] + D[48] +
338  D[47] + D[46] + D[45] + D[43] + D[41] + D[40] +
339  D[39] + D[38] + D[37] + D[36] + D[35] + D[34] +
340  D[33] + D[32] + D[31] + D[30] + D[27] + D[26] +
341  D[25] + D[24] + D[23] + D[22] + D[21] + D[20] +
342  D[19] + D[18] + D[17] + D[16] + D[15] + D[13] +
343  D[12] + D[11] + D[10] + D[9] + D[8] + D[7] +
344  D[6] + D[5] + D[4] + D[3] + D[2] + D[1] +
345  D[0] + C[0] + C[1] + C[2] + C[3] + C[4] +
346  C[5] + C[6] + C[7] + C[12] + C[13] + C[14] +
347  C[15] )%2;
348 
349  myCRC[1] = ( D[63] + D[62] + D[61] + D[56] + D[55] + D[54] +
350  D[53] + D[52] + D[51] + D[50] + D[49] + D[48] +
351  D[47] + D[46] + D[44] + D[42] + D[41] + D[40] +
352  D[39] + D[38] + D[37] + D[36] + D[35] + D[34] +
353  D[33] + D[32] + D[31] + D[28] + D[27] + D[26] +
354  D[25] + D[24] + D[23] + D[22] + D[21] + D[20] +
355  D[19] + D[18] + D[17] + D[16] + D[14] + D[13] +
356  D[12] + D[11] + D[10] + D[9] + D[8] + D[7] +
357  D[6] + D[5] + D[4] + D[3] + D[2] + D[1] +
358  C[0] + C[1] + C[2] + C[3] + C[4] + C[5] +
359  C[6] + C[7] + C[8] + C[13] + C[14] + C[15] )%2;
360 
361  myCRC[2] = ( D[61] + D[60] + D[57] + D[56] + D[46] + D[42] +
362  D[31] + D[30] + D[29] + D[28] + D[16] + D[14] +
363  D[1] + D[0] + C[8] + C[9] + C[12] + C[13] )%2;
364 
365  myCRC[3] = ( D[62] + D[61] + D[58] + D[57] + D[47] + D[43] +
366  D[32] + D[31] + D[30] + D[29] + D[17] + D[15] +
367  D[2] + D[1] + C[9] + C[10] + C[13] + C[14] )%2;
368 
369  myCRC[4] = ( D[63] + D[62] + D[59] + D[58] + D[48] + D[44] +
370  D[33] + D[32] + D[31] + D[30] + D[18] + D[16] +
371  D[3] + D[2] + C[0] + C[10] + C[11] + C[14] +
372  C[15] )%2;
373 
374  myCRC[5] = ( D[63] + D[60] + D[59] + D[49] + D[45] + D[34] +
375  D[33] + D[32] + D[31] + D[19] + D[17] + D[4] +
376  D[3] + C[1] + C[11] + C[12] + C[15] )%2;
377 
378  myCRC[6] = ( D[61] + D[60] + D[50] + D[46] + D[35] + D[34] +
379  D[33] + D[32] + D[20] + D[18] + D[5] + D[4] +
380  C[2] + C[12] + C[13] )%2;
381 
382  myCRC[7] = ( D[62] + D[61] + D[51] + D[47] + D[36] + D[35] +
383  D[34] + D[33] + D[21] + D[19] + D[6] + D[5] +
384  C[3] + C[13] + C[14] )%2;
385 
386  myCRC[8] = ( D[63] + D[62] + D[52] + D[48] + D[37] + D[36] +
387  D[35] + D[34] + D[22] + D[20] + D[7] + D[6] +
388  C[0] + C[4] + C[14] + C[15] )%2;
389 
390  myCRC[9] = ( D[63] + D[53] + D[49] + D[38] + D[37] + D[36] +
391  D[35] + D[23] + D[21] + D[8] + D[7] + C[1] +
392  C[5] + C[15] )%2;
393 
394  myCRC[10] = ( D[54] + D[50] + D[39] + D[38] + D[37] + D[36] +
395  D[24] + D[22] + D[9] + D[8] + C[2] + C[6] )%2;
396 
397  myCRC[11] = ( D[55] + D[51] + D[40] + D[39] + D[38] + D[37] +
398  D[25] + D[23] + D[10] + D[9] + C[3] + C[7] )%2;
399 
400  myCRC[12] = ( D[56] + D[52] + D[41] + D[40] + D[39] + D[38] +
401  D[26] + D[24] + D[11] + D[10] + C[4] + C[8] )%2;
402 
403  myCRC[13] = ( D[57] + D[53] + D[42] + D[41] + D[40] + D[39] +
404  D[27] + D[25] + D[12] + D[11] + C[5] + C[9] )%2;
405 
406  myCRC[14] = ( D[58] + D[54] + D[43] + D[42] + D[41] + D[40] +
407  D[28] + D[26] + D[13] + D[12] + C[6] + C[10] )%2;
408 
409  myCRC[15] = ( D[63] + D[62] + D[61] + D[60] + D[59] + D[54] +
410  D[53] + D[52] + D[51] + D[50] + D[49] + D[48] +
411  D[47] + D[46] + D[45] + D[44] + D[42] + D[40] +
412  D[39] + D[38] + D[37] + D[36] + D[35] + D[34] +
413  D[33] + D[32] + D[31] + D[30] + D[29] + D[26] +
414  D[25] + D[24] + D[23] + D[22] + D[21] + D[20] +
415  D[19] + D[18] + D[17] + D[16] + D[15] + D[14] +
416  D[12] + D[11] + D[10] + D[9] + D[8] + D[7] +
417  D[6] + D[5] + D[4] + D[3] + D[2] + D[1] +
418  D[0] + C[0] + C[1] + C[2] + C[3] + C[4] +
419  C[5] + C[6] + C[11] + C[12] + C[13] + C[14] +
420  C[15] )%2;
421 
422  int tempC = 0x0;
423  for(int i=0; i<16 ; i++){ tempC = tempC + (myCRC[i]<<i); }
424  myC = tempC ;
425  return;
426 }
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:44
int i
Definition: DBlmapReader.cc:9
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
virtual ~DTTFFEDSim()
Destructor.
Definition: DTTFFEDSim.cc:38
edm::EDGetTokenT< L1MuDTChambThContainer > ChTh_tok
Definition: DTTFFEDSim.h:68
The_Container::const_iterator The_iterator
int sector(int channel)
Definition: DTTFFEDSim.cc:304
unsigned int eventNum
Definition: DTTFFEDSim.h:48
const FEDRawData & FEDData(int fedid) const
retrieve data for fed
edm::InputTag DTDigiInputTag
Definition: DTTFFEDSim.h:50
void resize(size_t newsize)
Definition: FEDRawData.cc:32
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
int wheel(int channel)
Definition: DTTFFEDSim.cc:315
void calcCRC(int myD1, int myD2, int &myC)
Definition: DTTFFEDSim.cc:328
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Phi_Container::const_iterator Phi_iterator
int channel(int wheel, int sector, int bx)
Definition: DTTFFEDSim.cc:271
void produce(edm::Event &e, const edm::EventSetup &c)
Produce digis out of raw data.
Definition: DTTFFEDSim.cc:40
TrackContainer::const_iterator Trackiterator
DecomposeProduct< arg, typename Div::arg > D
Definition: Factorize.h:150
bool fillRawData(edm::Event &e, FEDRawDataCollection &data)
Generate and fill FED raw data for a full event.
Definition: DTTFFEDSim.cc:52
int bxNr(int channel)
Definition: DTTFFEDSim.cc:291
edm::EventID id() const
Definition: EventBase.h:56
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
edm::EDGetTokenT< L1MuDTTrackContainer > Trk_tok
Definition: DTTFFEDSim.h:69
edm::EDGetTokenT< L1MuDTChambPhContainer > ChPh_tok
Definition: DTTFFEDSim.h:67
const unsigned char * data() const
Return a const pointer to the beginning of the data buffer.
Definition: FEDRawData.cc:28
edm::InputTag DTPHTFInputTag
Definition: DTTFFEDSim.h:51
DTTFFEDSim(const edm::ParameterSet &pset)
Constructor.
Definition: DTTFFEDSim.cc:25