CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
DTTFFEDSim Class Reference

#include <DTTFFEDSim.h>

Inheritance diagram for DTTFFEDSim:
edm::stream::EDProducer<>

Public Member Functions

 DTTFFEDSim (const edm::ParameterSet &pset)
 Constructor. More...
 
bool fillRawData (edm::Event &e, FEDRawDataCollection &data)
 Generate and fill FED raw data for a full event. More...
 
void produce (edm::Event &e, const edm::EventSetup &c) override
 Produce digis out of raw data. More...
 
 ~DTTFFEDSim () override
 Destructor. More...
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Private Member Functions

int bxNr (int channel)
 
int channel (int wheel, int sector, int bx)
 
edm::InputTag getDTDigiInputTag ()
 
edm::InputTag getDTPHTFInputTag ()
 
int sector (int channel)
 
int wheel (int channel)
 

Private Attributes

edm::EDGetTokenT< L1MuDTChambPhContainerChPh_tok
 
edm::EDGetTokenT< L1MuDTChambThContainerChTh_tok
 
edm::InputTag DTDigiInputTag
 
edm::InputTag DTPHTFInputTag
 
unsigned int eventNum
 
edm::EDGetTokenT< L1MuDTTrackContainerTrk_tok
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

Detailed Description

L1 DT Track Finder Digi-to-Raw

J. Troconiz UAM Madrid

Definition at line 29 of file DTTFFEDSim.h.

Constructor & Destructor Documentation

◆ DTTFFEDSim()

DTTFFEDSim::DTTFFEDSim ( const edm::ParameterSet pset)

Constructor.

Definition at line 26 of file DTTFFEDSim.cc.

References ChPh_tok, ChTh_tok, DTDigiInputTag, DTPHTFInputTag, muonDTDigis_cfi::pset, and Trk_tok.

26  : eventNum(0) {
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 }
edm::EDGetTokenT< L1MuDTChambThContainer > ChTh_tok
Definition: DTTFFEDSim.h:62
unsigned int eventNum
Definition: DTTFFEDSim.h:44
edm::InputTag DTDigiInputTag
Definition: DTTFFEDSim.h:46
edm::EDGetTokenT< L1MuDTTrackContainer > Trk_tok
Definition: DTTFFEDSim.h:63
edm::EDGetTokenT< L1MuDTChambPhContainer > ChPh_tok
Definition: DTTFFEDSim.h:61
edm::InputTag DTPHTFInputTag
Definition: DTTFFEDSim.h:47

◆ ~DTTFFEDSim()

DTTFFEDSim::~DTTFFEDSim ( )
override

Destructor.

Definition at line 37 of file DTTFFEDSim.cc.

37 {}

Member Function Documentation

◆ bxNr()

int DTTFFEDSim::bxNr ( int  channel)
private

Definition at line 272 of file DTTFFEDSim.cc.

References channel().

272  {
273  int myChannel = channel;
274 
275  if (myChannel > 127)
276  myChannel -= 2;
277 
278  if (myChannel < 0 || myChannel > 251) {
279  return -999;
280  }
281 
282  int myBx = 1 - (myChannel % 3);
283 
284  return myBx;
285 }
int channel(int wheel, int sector, int bx)
Definition: DTTFFEDSim.cc:246

◆ channel()

int DTTFFEDSim::channel ( int  wheel,
int  sector,
int  bx 
)
private

Definition at line 246 of file DTTFFEDSim.cc.

References funct::abs(), nano_mu_digi_cff::bx, sector(), and wheel().

Referenced by bxNr(), fillRawData(), sector(), and wheel().

246  {
247  // wheel : -3 -2 -1 +1 +2 +3 <=> PHTF's : N2, N1, N0, P0, P1, P2
248  // 0 <=> ETTF
249  // sector : 0 -> 11
250  // bx : -1 -> +1
251 
252  int myChannel = 255;
253 
254  if (abs(bx) > 1) {
255  return myChannel;
256  }
257  if (sector < 0 || sector > 11) {
258  return myChannel;
259  }
260  if (abs(wheel) > 3) {
261  return myChannel;
262  }
263 
264  myChannel = sector * 21 + wheel * 3 - bx + 10;
265 
266  if (myChannel > 125)
267  myChannel += 2;
268 
269  return myChannel;
270 }
int sector(int channel)
Definition: DTTFFEDSim.cc:287
int wheel(int channel)
Definition: DTTFFEDSim.cc:300
Abs< T >::type abs(const T &t)
Definition: Abs.h:22

◆ fillRawData()

bool DTTFFEDSim::fillRawData ( edm::Event e,
FEDRawDataCollection data 
)

Generate and fill FED raw data for a full event.

Definition at line 50 of file DTTFFEDSim.cc.

References L1MuDTChambPhContainer::bxSize(), L1MuDTTrackContainer::bxSize(), L1MuDTChambThContainer::bxSize(), dt_crc::calcCRC(), channel(), ChPh_tok, ChTh_tok, FEDRawData::data(), data, MillePedeFileConverter_cfg::e, eventNum, L1MuDTChambPhContainer::getContainer(), L1MuDTTrackContainer::getContainer(), L1MuDTChambThContainer::getContainer(), createfilelist::int, groupFilesInBlocks::lines, FEDRawData::resize(), and Trk_tok.

Referenced by produce().

50  {
51  eventNum = e.id().event();
52 
53  int lines = 2;
54 
56  e.getByToken(ChPh_tok, phtrig);
57  lines += phtrig->bxSize(-1, 1);
58 
60  e.getByToken(ChTh_tok, thtrig);
61  lines += thtrig->bxSize(-1, 1);
62 
64  e.getByToken(Trk_tok, trtrig);
65  lines += trtrig->bxSize(-1, 1) * 3;
66 
67  FEDRawData &dttfdata = data.FEDData(0x30C);
68  dttfdata.resize(lines * 8); // size in bytes
69  unsigned char *LineFED = dttfdata.data();
70 
71  int *dataWord1 = new int;
72  int *dataWord2 = new int;
73 
74  //--> Header
75 
76  *dataWord1 = 0x50000000 + (eventNum & 0xFFFFFF);
77  *dataWord2 = 0x00030C00;
78 
79  int newCRC = 0xFFFF;
80  dt_crc::calcCRC(*dataWord1, *dataWord2, newCRC);
81 
82  *((int *)LineFED) = *dataWord2;
83  LineFED += 4;
84  *((int *)LineFED) = *dataWord1;
85 
86  //--> DTTF data
87 
88  int TS1Id[4], TS2Id[4]; // word identifier for TS #1,#2 for stations
89  TS1Id[0] = 0x0E;
90  TS2Id[0] = 0x1E;
91  TS1Id[1] = 0x2E;
92  TS2Id[1] = 0x3E;
93  TS1Id[3] = 0x4E;
94  TS2Id[3] = 0x5E;
95  TS1Id[2] = 0x8FFF8;
96  TS2Id[2] = 0x9FFF8;
97 
98  // Input
100 
101  for (tsphi = phtrig->getContainer()->begin(); tsphi != phtrig->getContainer()->end(); tsphi++) {
102  if (tsphi->code() != 7) {
103  int wheelID = tsphi->whNum() + 1;
104  if (wheelID <= 0)
105  wheelID -= 2;
106  int stationID = tsphi->stNum() - 1;
107  int is2nd = tsphi->Ts2Tag();
108 
109  int channelNr = channel(wheelID, tsphi->scNum(), tsphi->bxNum() - is2nd);
110  if (channelNr == 255)
111  continue;
112  int TSId = (is2nd == 0) ? TS1Id[stationID] : TS2Id[stationID];
113 
114  *dataWord1 = ((channelNr & 0xFF) << 24) + 0x00FFFFFF;
115 
116  if (stationID != 2) {
117  *dataWord2 = ((TSId & 0x0FF) << 24) + (~(tsphi->code() + 1) & 0x007) + ((~tsphi->phiB() & 0x3FF) << 3) +
118  ((~tsphi->phi() & 0xFFF) << 13);
119  } else {
120  *dataWord2 = ((TSId & 0xFFFFF) << 12) + (~(tsphi->code() + 1) & 0x00007) + ((~tsphi->phi() & 0x00FFF) << 3);
121  }
122 
123  dt_crc::calcCRC(*dataWord1, *dataWord2, newCRC);
124 
125  LineFED += 4;
126  *((int *)LineFED) = *dataWord2;
127  LineFED += 4;
128  *((int *)LineFED) = *dataWord1;
129  }
130  }
131  // Input
132 
133  // Input
135 
136  for (tsthe = thtrig->getContainer()->begin(); tsthe != thtrig->getContainer()->end(); tsthe++) {
137  int wheelTh = tsthe->whNum();
138  int sectorID = tsthe->scNum();
139 
140  int channelNr = channel(0, sectorID, tsthe->bxNum());
141  if (channelNr == 255)
142  continue;
143  int TSId = wheelTh + 2;
144 
145  *dataWord1 = ((channelNr & 0xFF) << 24) + 0x00FFFFFF;
146 
147  *dataWord2 = ((TSId & 0x07) << 28) + 0x0FFFFFFF;
148 
149  int stationID = tsthe->stNum() - 1;
150  for (int bti = 0; bti < 7; bti++)
151  if (wheelTh == -2 || wheelTh == -1 ||
152  (wheelTh == 0 &&
153  (sectorID == 0 || sectorID == 3 || sectorID == 4 || sectorID == 7 || sectorID == 8 || sectorID == 11)))
154  *dataWord2 -= (tsthe->position(bti) & 0x1) << (stationID * 7 + bti);
155  else
156  *dataWord2 -= (tsthe->position(6 - bti) & 0x1) << (stationID * 7 + bti);
157 
158  dt_crc::calcCRC(*dataWord1, *dataWord2, newCRC);
159 
160  LineFED += 4;
161  *((int *)LineFED) = *dataWord2;
162  LineFED += 4;
163  *((int *)LineFED) = *dataWord1;
164  }
165  // Input
166 
167  // Output
169 
170  for (tstrk = trtrig->getContainer()->begin(); tstrk != trtrig->getContainer()->end(); tstrk++) {
171  int channelNr = channel(tstrk->whNum(), tstrk->scNum(), tstrk->bx());
172  if (channelNr == 255)
173  continue;
174  int TSId = (tstrk->TrkTag() == 0) ? 0xAFFF : 0xBFFF;
175 
176  *dataWord1 = ((channelNr & 0xFF) << 24) + 0x00FFFFFF;
177 
178  *dataWord2 = ((TSId & 0xFFFF) << 16) + (tstrk->stNum(4) & 0x0000F) + ((tstrk->stNum(3) & 0x0000F) << 4) +
179  ((tstrk->stNum(2) & 0x0000F) << 8) + ((tstrk->stNum(1) & 0x00003) << 12);
180 
181  dt_crc::calcCRC(*dataWord1, *dataWord2, newCRC);
182 
183  LineFED += 4;
184  *((int *)LineFED) = *dataWord2;
185  LineFED += 4;
186  *((int *)LineFED) = *dataWord1;
187 
188  TSId = (tstrk->TrkTag() == 0) ? 0xCFFE : 0xDFFE;
189 
190  *dataWord1 = ((channelNr & 0xFF) << 24) + 0x00FFFFFF;
191 
192  *dataWord2 = ((TSId & 0xFFFE) << 16) + (~tstrk->quality_packed() & 0x0007) + ((tstrk->phi_packed() & 0x00FF) << 3) +
193  ((~tstrk->charge_packed() & 0x0001) << 11) + ((~tstrk->pt_packed() & 0x001F) << 12);
194 
195  dt_crc::calcCRC(*dataWord1, *dataWord2, newCRC);
196 
197  LineFED += 4;
198  *((int *)LineFED) = *dataWord2;
199  LineFED += 4;
200  *((int *)LineFED) = *dataWord1;
201 
202  channelNr = channel(0, tstrk->scNum(), tstrk->bx());
203  if (channelNr == 255)
204  continue;
205  TSId = (tstrk->whNum() + 3) << 16;
206  TSId += (tstrk->whNum() < 0) ? 0x8FFFC : 0x7FFFC;
207 
208  *dataWord1 = ((channelNr & 0xFF) << 24) + 0x00FFFFFF;
209 
210  *dataWord2 = (TSId & 0xFFFFC) << 12;
211 
212  if (tstrk->TrkTag() == 0) {
213  *dataWord2 += 0x3F80 + (tstrk->eta_packed() & 0x003F) + ((~tstrk->finehalo_packed() & 0x0001) << 6);
214  } else {
215  *dataWord2 += 0x007F + ((tstrk->eta_packed() & 0x003F) << 7) + ((~tstrk->finehalo_packed() & 0x0001) << 13);
216  }
217 
218  dt_crc::calcCRC(*dataWord1, *dataWord2, newCRC);
219 
220  LineFED += 4;
221  *((int *)LineFED) = *dataWord2;
222  LineFED += 4;
223  *((int *)LineFED) = *dataWord1;
224  }
225  // Output
226 
227  //--> Trailer
228 
229  *dataWord1 = 0xA0000000 + (lines & 0xFFFFFF);
230  *dataWord2 = 0;
231 
232  dt_crc::calcCRC(*dataWord1, *dataWord2 & 0xFFFF, newCRC);
233 
234  *dataWord2 += (newCRC & 0xFFFF) << 16;
235 
236  LineFED += 4;
237  *((int *)LineFED) = *dataWord2;
238  LineFED += 4;
239  *((int *)LineFED) = *dataWord1;
240 
241  delete dataWord1;
242  delete dataWord2;
243  return true;
244 }
int bxSize(int step1, int step2) const
int bxSize(int step1, int step2) const
edm::EDGetTokenT< L1MuDTChambThContainer > ChTh_tok
Definition: DTTFFEDSim.h:62
The_Container::const_iterator The_iterator
unsigned int eventNum
Definition: DTTFFEDSim.h:44
The_Container const * getContainer() const
void resize(size_t newsize)
Definition: FEDRawData.cc:28
TrackContainer const * getContainer() const
Phi_Container::const_iterator Phi_iterator
int channel(int wheel, int sector, int bx)
Definition: DTTFFEDSim.cc:246
TrackContainer::const_iterator Trackiterator
Phi_Container const * getContainer() const
void calcCRC(long, int &)
Definition: DTCRC.cc:3
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:80
edm::EDGetTokenT< L1MuDTTrackContainer > Trk_tok
Definition: DTTFFEDSim.h:63
int bxSize(int step1, int step2) const
edm::EDGetTokenT< L1MuDTChambPhContainer > ChPh_tok
Definition: DTTFFEDSim.h:61
const unsigned char * data() const
Return a const pointer to the beginning of the data buffer.
Definition: FEDRawData.cc:24

◆ getDTDigiInputTag()

edm::InputTag DTTFFEDSim::getDTDigiInputTag ( )
inlineprivate

Definition at line 58 of file DTTFFEDSim.h.

References DTDigiInputTag.

58 { return DTDigiInputTag; }
edm::InputTag DTDigiInputTag
Definition: DTTFFEDSim.h:46

◆ getDTPHTFInputTag()

edm::InputTag DTTFFEDSim::getDTPHTFInputTag ( )
inlineprivate

Definition at line 59 of file DTTFFEDSim.h.

References DTPHTFInputTag.

59 { return DTPHTFInputTag; }
edm::InputTag DTPHTFInputTag
Definition: DTTFFEDSim.h:47

◆ produce()

void DTTFFEDSim::produce ( edm::Event e,
const edm::EventSetup c 
)
override

Produce digis out of raw data.

Definition at line 39 of file DTTFFEDSim.cc.

References data, MillePedeFileConverter_cfg::e, EcalFEDMonitor_cfi::FEDRawDataCollection, fillRawData(), and eostools::move().

39  {
41 
42  if (!fillRawData(e, data))
43  return;
44 
45  unique_ptr<FEDRawDataCollection> fed_product(new FEDRawDataCollection(data));
46 
47  e.put(std::move(fed_product));
48 }
bool fillRawData(edm::Event &e, FEDRawDataCollection &data)
Generate and fill FED raw data for a full event.
Definition: DTTFFEDSim.cc:50
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:80
def move(src, dest)
Definition: eostools.py:511

◆ sector()

int DTTFFEDSim::sector ( int  channel)
private

Definition at line 287 of file DTTFFEDSim.cc.

References channel().

Referenced by channel(), and geometryXMLparser.DTAlignable::index().

287  {
288  int myChannel = channel;
289 
290  if (myChannel > 127)
291  myChannel -= 2;
292 
293  if (myChannel < 0 || myChannel > 251) {
294  return -999;
295  }
296 
297  return myChannel / 21;
298 }
int channel(int wheel, int sector, int bx)
Definition: DTTFFEDSim.cc:246

◆ wheel()

int DTTFFEDSim::wheel ( int  channel)
private

Definition at line 300 of file DTTFFEDSim.cc.

References channel().

Referenced by channel(), and geometryXMLparser.DTAlignable::index().

300  {
301  int myChannel = channel;
302 
303  if (myChannel > 127)
304  myChannel -= 2;
305 
306  if (myChannel < 0 || myChannel > 251) {
307  return -999;
308  }
309 
310  int myWheel = ((myChannel % 21) / 3) - 3;
311 
312  return myWheel;
313 }
int channel(int wheel, int sector, int bx)
Definition: DTTFFEDSim.cc:246

Member Data Documentation

◆ ChPh_tok

edm::EDGetTokenT<L1MuDTChambPhContainer> DTTFFEDSim::ChPh_tok
private

Definition at line 61 of file DTTFFEDSim.h.

Referenced by DTTFFEDSim(), and fillRawData().

◆ ChTh_tok

edm::EDGetTokenT<L1MuDTChambThContainer> DTTFFEDSim::ChTh_tok
private

Definition at line 62 of file DTTFFEDSim.h.

Referenced by DTTFFEDSim(), and fillRawData().

◆ DTDigiInputTag

edm::InputTag DTTFFEDSim::DTDigiInputTag
private

Definition at line 46 of file DTTFFEDSim.h.

Referenced by DTTFFEDSim(), and getDTDigiInputTag().

◆ DTPHTFInputTag

edm::InputTag DTTFFEDSim::DTPHTFInputTag
private

Definition at line 47 of file DTTFFEDSim.h.

Referenced by DTTFFEDSim(), and getDTPHTFInputTag().

◆ eventNum

unsigned int DTTFFEDSim::eventNum
private

Definition at line 44 of file DTTFFEDSim.h.

Referenced by fillRawData().

◆ Trk_tok

edm::EDGetTokenT<L1MuDTTrackContainer> DTTFFEDSim::Trk_tok
private

Definition at line 63 of file DTTFFEDSim.h.

Referenced by DTTFFEDSim(), and fillRawData().