CMS 3D CMS Logo

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 
21 
22 #include <iostream>
23 
24 using namespace std;
25 
27 
28  produces<FEDRawDataCollection>();
29 
30  DTDigiInputTag = pset.getParameter<edm::InputTag>("DTDigi_Source");
31  DTPHTFInputTag = pset.getParameter<edm::InputTag>("DTTracks_Source");
32 
33  ChPh_tok = consumes<L1MuDTChambPhContainer>(DTDigiInputTag);
34  ChTh_tok = consumes<L1MuDTChambThContainer>(DTDigiInputTag);
35  Trk_tok = consumes<L1MuDTTrackContainer>(DTPHTFInputTag);
36 
37 }
38 
40 
42 
44 
45  if (!fillRawData(e, data)) return;
46 
47  unique_ptr<FEDRawDataCollection> fed_product(new FEDRawDataCollection(data));
48 
49  e.put(std::move(fed_product));
50 
51 }
52 
55 
56  eventNum = e.id().event();
57 
58  int lines = 2;
59 
61  e.getByToken(ChPh_tok,phtrig);
62  lines += phtrig->bxSize(-1, 1);
63 
65  e.getByToken(ChTh_tok,thtrig);
66  lines += thtrig->bxSize(-1, 1);
67 
69  e.getByToken(Trk_tok,trtrig);
70  lines += trtrig->bxSize(-1, 1)*3;
71 
72  FEDRawData& dttfdata = data.FEDData(0x30C);
73  dttfdata.resize(lines*8); // size in bytes
74  unsigned char* LineFED=dttfdata.data();
75 
76  int* dataWord1 = new int;
77  int* dataWord2 = new int;
78 
79  //--> Header
80 
81  *dataWord1 = 0x50000000
82  + (eventNum&0xFFFFFF);
83  *dataWord2 = 0x00030C00;
84 
85  int newCRC = 0xFFFF;
86  dt_crc::calcCRC(*dataWord1, *dataWord2, newCRC);
87 
88  *((int*)LineFED)=*dataWord2;
89  LineFED+=4;
90  *((int*)LineFED)=*dataWord1;
91 
92  //--> DTTF data
93 
94  int TS1Id[4], TS2Id[4]; // word identifier for TS #1,#2 for stations
95  TS1Id[0] = 0x0E;
96  TS2Id[0] = 0x1E;
97  TS1Id[1] = 0x2E;
98  TS2Id[1] = 0x3E;
99  TS1Id[3] = 0x4E;
100  TS2Id[3] = 0x5E;
101  TS1Id[2] = 0x8FFF8;
102  TS2Id[2] = 0x9FFF8;
103 
104  //Input
106 
107  for ( tsphi = phtrig->getContainer()->begin();
108  tsphi != phtrig->getContainer()->end();
109  tsphi++ ) {
110  if ( tsphi->code() != 7 ) {
111 
112  int wheelID = tsphi->whNum()+1;
113  if ( wheelID <= 0 ) wheelID -= 2;
114  int stationID = tsphi->stNum()-1;
115  int is2nd = tsphi->Ts2Tag();
116 
117  int channelNr = channel(wheelID, tsphi->scNum(), tsphi->bxNum()-is2nd);
118  if ( channelNr == 255 ) continue;
119  int TSId = ( is2nd == 0 ) ? TS1Id[stationID] : TS2Id[stationID];
120 
121  *dataWord1 = ((channelNr&0xFF)<<24)
122  + 0x00FFFFFF;
123 
124  if ( stationID != 2 ){
125  *dataWord2 = ( (TSId&0x0FF)<<24)
126  + (~(tsphi->code()+1)&0x007)
127  + ( (~tsphi->phiB()&0x3FF)<<3)
128  + ( (~tsphi->phi()&0xFFF)<<13);
129  }
130  else {
131  *dataWord2 = ( (TSId&0xFFFFF)<<12)
132  + (~(tsphi->code()+1)&0x00007)
133  + ( (~tsphi->phi()&0x00FFF)<<3);
134  }
135 
136  dt_crc::calcCRC(*dataWord1, *dataWord2, newCRC);
137 
138  LineFED+=4;
139  *((int*)LineFED)=*dataWord2;
140  LineFED+=4;
141  *((int*)LineFED)=*dataWord1;
142  }
143  }
144  //Input
145 
146  //Input
148 
149  for ( tsthe = thtrig->getContainer()->begin();
150  tsthe != thtrig->getContainer()->end();
151  tsthe++ ) {
152 
153  int wheelTh = tsthe->whNum();
154  int sectorID = tsthe->scNum();
155 
156  int channelNr = channel(0, sectorID, tsthe->bxNum());
157  if ( channelNr == 255 ) continue;
158  int TSId = wheelTh+2;
159 
160  *dataWord1 = ((channelNr&0xFF)<<24)
161  + 0x00FFFFFF;
162 
163  *dataWord2 = ((TSId&0x07)<<28)
164  + 0x0FFFFFFF;
165 
166  int stationID = tsthe->stNum()-1;
167  for ( int bti = 0; bti < 7; bti++ )
168  if ( wheelTh == -2 || wheelTh == -1 ||
169  ( wheelTh == 0 && (sectorID == 0 || sectorID == 3 || sectorID == 4 || sectorID == 7 || sectorID == 8 || sectorID == 11) ) )
170  *dataWord2 -= (tsthe->position(bti)&0x1)<<(stationID*7+bti);
171  else
172  *dataWord2 -= (tsthe->position(6-bti)&0x1)<<(stationID*7+bti);
173 
174  dt_crc::calcCRC(*dataWord1, *dataWord2, newCRC);
175 
176  LineFED+=4;
177  *((int*)LineFED)=*dataWord2;
178  LineFED+=4;
179  *((int*)LineFED)=*dataWord1;
180  }
181  //Input
182 
183  //Output
185 
186  for ( tstrk = trtrig->getContainer()->begin();
187  tstrk != trtrig->getContainer()->end();
188  tstrk++ ) {
189 
190  int channelNr = channel(tstrk->whNum(), tstrk->scNum(), tstrk->bx());
191  if ( channelNr == 255 ) continue;
192  int TSId = ( tstrk->TrkTag() == 0 ) ? 0xAFFF : 0xBFFF;
193 
194  *dataWord1 = ((channelNr&0xFF)<<24)
195  + 0x00FFFFFF;
196 
197  *dataWord2 = ( (TSId&0xFFFF)<<16)
198  + ( tstrk->stNum(4)&0x0000F)
199  + ((tstrk->stNum(3)&0x0000F)<<4)
200  + ((tstrk->stNum(2)&0x0000F)<<8)
201  + ((tstrk->stNum(1)&0x00003)<<12);
202 
203  dt_crc::calcCRC(*dataWord1, *dataWord2, newCRC);
204 
205  LineFED+=4;
206  *((int*)LineFED)=*dataWord2;
207  LineFED+=4;
208  *((int*)LineFED)=*dataWord1;
209 
210  TSId = ( tstrk->TrkTag() == 0 ) ? 0xCFFE : 0xDFFE;
211 
212  *dataWord1 = ((channelNr&0xFF)<<24)
213  + 0x00FFFFFF;
214 
215  *dataWord2 = ( (TSId&0xFFFE)<<16)
216  + ( ~tstrk->quality_packed()&0x0007)
217  + ( (tstrk->phi_packed()&0x00FF)<<3)
218  + ( (~tstrk->charge_packed()&0x0001)<<11)
219  + ( (~tstrk->pt_packed()&0x001F)<<12);
220 
221  dt_crc::calcCRC(*dataWord1, *dataWord2, newCRC);
222 
223  LineFED+=4;
224  *((int*)LineFED)=*dataWord2;
225  LineFED+=4;
226  *((int*)LineFED)=*dataWord1;
227 
228  channelNr = channel(0, tstrk->scNum(), tstrk->bx());
229  if ( channelNr == 255 ) continue;
230  TSId = (tstrk->whNum()+3)<<16;
231  TSId += ( tstrk->whNum() < 0 ) ? 0x8FFFC : 0x7FFFC;
232 
233  *dataWord1 = ((channelNr&0xFF)<<24)
234  + 0x00FFFFFF;
235 
236  *dataWord2 = (TSId&0xFFFFC)<<12;
237 
238  if ( tstrk->TrkTag() == 0 ) {
239  *dataWord2 += 0x3F80
240  + ( tstrk->eta_packed()&0x003F)
241  + ((~tstrk->finehalo_packed()&0x0001)<<6);
242  }
243  else {
244  *dataWord2 += 0x007F
245  + ( ( tstrk->eta_packed()&0x003F)<<7)
246  + ((~tstrk->finehalo_packed()&0x0001)<<13);
247  }
248 
249  dt_crc::calcCRC(*dataWord1, *dataWord2, newCRC);
250 
251  LineFED+=4;
252  *((int*)LineFED)=*dataWord2;
253  LineFED+=4;
254  *((int*)LineFED)=*dataWord1;
255  }
256  //Output
257 
258  //--> Trailer
259 
260  *dataWord1 = 0xA0000000
261  + (lines&0xFFFFFF);
262  *dataWord2 = 0;
263 
264  dt_crc::calcCRC(*dataWord1, *dataWord2&0xFFFF, newCRC);
265 
266  *dataWord2 += (newCRC&0xFFFF)<<16;
267 
268  LineFED+=4;
269  *((int*)LineFED)=*dataWord2;
270  LineFED+=4;
271  *((int*)LineFED)=*dataWord1;
272 
273  delete dataWord1;
274  delete dataWord2;
275  return true;
276 }
277 
278 int DTTFFEDSim::channel( int wheel, int sector, int bx ){
279 
280  // wheel : -3 -2 -1 +1 +2 +3 <=> PHTF's : N2, N1, N0, P0, P1, P2
281  // 0 <=> ETTF
282  // sector : 0 -> 11
283  // bx : -1 -> +1
284 
285  int myChannel = 255;
286 
287  if ( abs(bx) > 1) { return myChannel; }
288  if ( sector < 0 || sector > 11) { return myChannel; }
289  if ( abs(wheel) > 3) { return myChannel; }
290 
291  myChannel = sector*21 + wheel*3 - bx + 10 ;
292 
293  if (myChannel > 125) myChannel += 2;
294 
295  return myChannel;
296 }
297 
299 
300  int myChannel = channel;
301 
302  if (myChannel > 127) myChannel -= 2;
303 
304  if (myChannel < 0 || myChannel > 251 ){ return -999; }
305 
306  int myBx = 1-(myChannel%3);
307 
308  return myBx;
309 }
310 
312 
313  int myChannel = channel;
314 
315  if (myChannel > 127) myChannel -= 2;
316 
317  if (myChannel < 0 || myChannel > 251 ){ return -999; }
318 
319  return myChannel/21;
320 }
321 
323 
324  int myChannel = channel;
325 
326  if (myChannel > 127) myChannel -= 2;
327 
328  if (myChannel < 0 || myChannel > 251 ){ return -999; }
329 
330  int myWheel = ((myChannel%21)/3)-3;
331 
332  return myWheel;
333 }
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:41
The_Container const * getContainer() const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
void produce(edm::Event &e, const edm::EventSetup &c) override
Produce digis out of raw data.
Definition: DTTFFEDSim.cc:41
int bxSize(int step1, int step2) const
edm::EDGetTokenT< L1MuDTChambThContainer > ChTh_tok
Definition: DTTFFEDSim.h:66
The_Container::const_iterator The_iterator
int bxSize(int step1, int step2) const
int sector(int channel)
Definition: DTTFFEDSim.cc:311
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
int wheel(int channel)
Definition: DTTFFEDSim.cc:322
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:278
TrackContainer::const_iterator Trackiterator
bool fillRawData(edm::Event &e, FEDRawDataCollection &data)
Generate and fill FED raw data for a full event.
Definition: DTTFFEDSim.cc:53
Phi_Container const * getContainer() const
int bxSize(int step1, int step2) const
~DTTFFEDSim() override
Destructor.
Definition: DTTFFEDSim.cc:39
int bxNr(int channel)
Definition: DTTFFEDSim.cc:298
edm::EventID id() const
Definition: EventBase.h:60
void calcCRC(long, int &)
Definition: DTCRC.cc:3
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
edm::EDGetTokenT< L1MuDTTrackContainer > Trk_tok
Definition: DTTFFEDSim.h:67
edm::EDGetTokenT< L1MuDTChambPhContainer > ChPh_tok
Definition: DTTFFEDSim.h:65
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
TrackContainer const * getContainer() const
def move(src, dest)
Definition: eostools.py:510
DTTFFEDSim(const edm::ParameterSet &pset)
Constructor.
Definition: DTTFFEDSim.cc:26